Galerkin FEM for elliptic PDEs
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.

FEM4EllipticPDE.cpp 37KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934
  1. // ===========================================================================
  2. //
  3. // Filename: FEM4EllipticPDE.cpp
  4. //
  5. // Created: 08/16/12 18:19:57
  6. // Compiler: Tested with g++, icpc, and MSVC 2010
  7. //
  8. // Author: Trevor Irons (ti)
  9. //
  10. // Organisation: Colorado School of Mines (CSM)
  11. // United States Geological Survey (USGS)
  12. //
  13. // Email: tirons@mines.edu, tirons@usgs.gov
  14. //
  15. // This program is free software: you can redistribute it and/or modify
  16. // it under the terms of the GNU General Public License as published by
  17. // the Free Software Foundation, either version 3 of the License, or
  18. // (at your option) any later version.
  19. //
  20. // This program is distributed in the hope that it will be useful,
  21. // but WITHOUT ANY WARRANTY; without even the implied warranty of
  22. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  23. // GNU General Public License for more details.
  24. //
  25. // You should have received a copy of the GNU General Public License
  26. // along with this program. If not, see <http://www.gnu.org/licenses/>.
  27. //
  28. // ===========================================================================
  29. /**
  30. @file
  31. @author Trevor Irons
  32. @date 08/16/12
  33. @version 0.0
  34. **/
  35. #include "FEM4EllipticPDE.h"
  36. namespace Lemma {
  37. std::ostream &operator<<(std::ostream &stream,
  38. const FEM4EllipticPDE &ob) {
  39. stream << *(LemmaObject*)(&ob);
  40. return stream;
  41. }
  42. // ==================== LIFECYCLE =======================
  43. FEM4EllipticPDE::FEM4EllipticPDE(const std::string&name) :
  44. LemmaObject(name), BndryH(1), BndrySigma(1),
  45. vtkSigma(NULL), vtkG(NULL), vtkGrid(NULL), gFcn3(NULL) {
  46. }
  47. FEM4EllipticPDE::~FEM4EllipticPDE() {
  48. }
  49. void FEM4EllipticPDE::Release() {
  50. delete this;
  51. }
  52. FEM4EllipticPDE* FEM4EllipticPDE::New( ) {
  53. FEM4EllipticPDE* Obj = new FEM4EllipticPDE("FEM4EllipticPDE");
  54. Obj->AttachTo(Obj);
  55. return Obj;
  56. }
  57. void FEM4EllipticPDE::Delete() {
  58. this->DetachFrom(this);
  59. }
  60. // ==================== OPERATIONS =======================
  61. void FEM4EllipticPDE::SetSigmaFunction(vtkImplicitFunction* sigma) {
  62. vtkSigma = sigma;
  63. }
  64. void FEM4EllipticPDE::SetBoundaryStep(const Real& h) {
  65. BndryH = h;
  66. }
  67. void FEM4EllipticPDE::SetGFunction(vtkImplicitFunction* g) {
  68. vtkG = g;
  69. }
  70. void FEM4EllipticPDE::SetGFunction( Real (*gFcn)(const Real&, const Real&, const Real&) ) {
  71. // vtkG = g;
  72. gFcn3 = gFcn;
  73. }
  74. void FEM4EllipticPDE::SetGrid(vtkDataSet* grid) {
  75. vtkGrid = grid;
  76. }
  77. vtkSmartPointer<vtkIdList> FEM4EllipticPDE::GetConnectedPoints(const int& id0) {
  78. vtkSmartPointer<vtkIdList> pointIds = vtkSmartPointer<vtkIdList>::New();
  79. vtkSmartPointer<vtkIdList> cellList = vtkSmartPointer<vtkIdList>::New();
  80. vtkGrid->GetPointCells(id0, cellList);
  81. for(int i=0;i<cellList->GetNumberOfIds(); ++i){
  82. vtkCell* cell = vtkGrid->GetCell(cellList->GetId(i));
  83. if(cell->GetNumberOfEdges() > 0){
  84. for(int j=0; j<cell->GetNumberOfEdges(); ++j){
  85. vtkCell* edge = cell->GetEdge(j);
  86. vtkIdList* edgePoints=edge->GetPointIds();
  87. if(edgePoints->GetId(0)==id0){
  88. pointIds->InsertUniqueId(edgePoints->GetId(1));
  89. } else if(edgePoints->GetId(1)==id0){
  90. pointIds->InsertUniqueId(edgePoints->GetId(0));
  91. }
  92. }
  93. }
  94. }
  95. return pointIds;
  96. }
  97. Real FEM4EllipticPDE::dist(Real r0[3], Real r1[3]) {
  98. Real rm0 = r1[0] - r0[0];
  99. Real rm1 = r1[1] - r0[1];
  100. Real rm2 = r1[2] - r0[2];
  101. return std::sqrt( rm0*rm0 + rm1*rm1 + rm2*rm2 );
  102. }
  103. Real FEM4EllipticPDE::dist(const Vector3r& r0, const Vector3r& r1) {
  104. Real rm0 = r1[0] - r0[0];
  105. Real rm1 = r1[1] - r0[1];
  106. Real rm2 = r1[2] - r0[2];
  107. return std::sqrt( rm0*rm0 + rm1*rm1 + rm2*rm2 );
  108. }
  109. //--------------------------------------------------------------------------------------
  110. // Class: FEM4EllipticPDE
  111. // Method: SetupDC
  112. //--------------------------------------------------------------------------------------
  113. void FEM4EllipticPDE::SetupDC ( DCSurvey* Survey, const int& ij ) {
  114. ////////////////////////////////////////////////////////////
  115. // Load vector g, solution vector u
  116. std::cout << "\nBuilding load vector (g)" << std::endl;
  117. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  118. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  119. int iia(0);
  120. Real jja(0);
  121. Survey->GetA( ij, iia, jja );
  122. //g(ii) = jj;
  123. int iib(0);
  124. Real jjb(0);
  125. Survey->GetB( ij, iib, jjb );
  126. //g(ii) = jj;
  127. /* 3D Phi */
  128. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  129. // Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  130. // for (int ip=0; ip<4; ++ip) {
  131. // double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  132. // C(ip, 0) = 1;
  133. // C(ip, 1) = pts[0];
  134. // C(ip, 2) = pts[1];
  135. // C(ip, 3) = pts[2];
  136. // }
  137. // Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  138. //
  139. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  140. int ID[4];
  141. ID[0] = Ids->GetId(0);
  142. ID[1] = Ids->GetId(1);
  143. ID[2] = Ids->GetId(2);
  144. ID[3] = Ids->GetId(3);
  145. //Real V = C.determinant(); // volume of tetrahedra
  146. Real sum(0);
  147. if (ID[0] == iia || ID[1] == iia || ID[2] == iia || ID[3] == iia ) {
  148. std::cout << "Caught A electrode, injecting " << iia << std::endl;
  149. //sum = 10;
  150. //g(ID[iia]) += jja;
  151. g(iia) += jja;
  152. }
  153. if (ID[0] == iib || ID[1] == iib || ID[2] == iib || ID[3] == iib) {
  154. //sum = -10;
  155. std::cout << "Caught B electrode, injecting " << iib << std::endl;
  156. //g(ID[iib]) += jjb;
  157. g(iib) += jjb;
  158. }
  159. //g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  160. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  161. }
  162. return ;
  163. } // ----- end of method FEM4EllipticPDE::SetupDC -----
  164. void FEM4EllipticPDE::Solve( const std::string& resfile) {
  165. ConstructAMatrix();
  166. //ConstructLoadVector();
  167. std::cout << "\nSolving" << std::endl;
  168. ////////////////////////////////////////////////////////////
  169. // Solving:
  170. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  171. //VectorXr u = chol.solve(g);
  172. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  173. Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  174. //Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  175. cg.setMaxIterations(3000);
  176. std::cout << "A: " << A.rows() << "\t" << A.cols() << std::endl;
  177. std::cout << "g: " << g.rows() << "\t" << g.cols() << std::endl;
  178. VectorXr u = cg.solve(g);
  179. std::cout << "#iterations: " << cg.iterations() << std::endl;
  180. std::cout << "estimated error: " << cg.error() << std::endl;
  181. vtkDoubleArray *gArray = vtkDoubleArray::New();
  182. vtkDoubleArray *uArray = vtkDoubleArray::New();
  183. uArray->SetNumberOfComponents(1);
  184. gArray->SetNumberOfComponents(1);
  185. for (int iu = 0; iu<u.size(); ++iu) {
  186. uArray->InsertTuple1(iu, u[iu]);
  187. gArray->InsertTuple1(iu, g[iu]);
  188. }
  189. uArray->SetName("u");
  190. gArray->SetName("g");
  191. vtkGrid->GetPointData()->AddArray(uArray);
  192. vtkGrid->GetPointData()->AddArray(gArray);
  193. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  194. Writer->SetInputData(vtkGrid);
  195. Writer->SetFileName(resfile.c_str());
  196. Writer->Write();
  197. Writer->Delete();
  198. gArray->Delete();
  199. uArray->Delete();
  200. }
  201. //--------------------------------------------------------------------------------------
  202. // Class: FEM4EllipticPDE
  203. // Method: ConstructAMatrix
  204. //--------------------------------------------------------------------------------------
  205. void FEM4EllipticPDE::ConstructAMatrix ( ) {
  206. /////////////////////////////////////////////////////////////////////////
  207. // Build stiffness matrix (A)
  208. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  209. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  210. //Eigen::SparseMatrix<Real>
  211. A.resize(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  212. std::vector< Eigen::Triplet<Real> > coeffs;
  213. if ( !vtkGrid->GetPointData()->GetScalars("vtkValidPointMask") ) {
  214. throw std::runtime_error("No vtkValidPointMask");
  215. }
  216. if ( !vtkGrid->GetCellData()->GetScalars("G") && !vtkGrid->GetPointData()->GetScalars("G") ) {
  217. throw std::runtime_error("No Cell or Point Data G");
  218. }
  219. bool GCell = false;
  220. if ( vtkGrid->GetCellData()->GetScalars("G") ) {
  221. GCell = true;
  222. }
  223. // Here we iterate over all of the cells
  224. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  225. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  226. // TODO, in production code we might not want to do this check here
  227. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  228. throw std::runtime_error("Non-tetrahedral mesh encountered!");
  229. }
  230. // construct coordinate matrix C
  231. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  232. for (int ip=0; ip<4; ++ip) {
  233. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  234. C(ip, 0) = 1;
  235. C(ip, 1) = pts[0] ;
  236. C(ip, 2) = pts[1] ;
  237. C(ip, 3) = pts[2] ;
  238. }
  239. Eigen::Matrix<Real, 4, 4> Phi = C.inverse(); // nabla \phi
  240. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  241. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  242. int ID[4];
  243. ID[0] = Ids->GetId(0);
  244. ID[1] = Ids->GetId(1);
  245. ID[2] = Ids->GetId(2);
  246. ID[3] = Ids->GetId(3);
  247. Real sum(0), sigma_bar(0);
  248. if (GCell) {
  249. sigma_bar = vtkGrid->GetCellData()->GetScalars("G")->GetTuple1(ic);
  250. } else {
  251. sigma_bar = vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0]);
  252. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1]);
  253. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2]);
  254. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3]);
  255. sigma_bar /= 4.;
  256. }
  257. for (int ip=0; ip<4; ++ip) {
  258. for (int ip2=0; ip2<4; ++ip2) {
  259. if (ip2 == ip) {
  260. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  261. // solve for the boundaries? Is one better? This seems to work, which is nice.
  262. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ip] ); // + sum;
  263. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ip])[0];
  264. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  265. coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], bdry + Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  266. } else {
  267. coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  268. }
  269. // Stiffness matrix no longer contains boundary conditions...
  270. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  271. }
  272. }
  273. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  274. }
  275. A.setFromTriplets(coeffs.begin(), coeffs.end());
  276. }
  277. void FEM4EllipticPDE::SetupPotential() {
  278. ////////////////////////////////////////////////////////////
  279. // Load vector g
  280. std::cout << "\nBuilding load vector (g)" << std::endl;
  281. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  282. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " cells" << std::endl;
  283. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  284. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  285. for (int ip=0; ip<4; ++ip) {
  286. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  287. C(ip, 0) = 1;
  288. C(ip, 1) = pts[0];
  289. C(ip, 2) = pts[1];
  290. C(ip, 3) = pts[2];
  291. }
  292. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  293. //Real V = C.determinant(); // volume of tetrahedra
  294. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  295. int ID[4];
  296. ID[0] = Ids->GetId(0);
  297. ID[1] = Ids->GetId(1);
  298. ID[2] = Ids->GetId(2);
  299. ID[3] = Ids->GetId(3);
  300. for (int ip=0; ip<4; ++ip) {
  301. g(ID[ip]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0] ) ;
  302. //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0]) > 1e-3 )
  303. }
  304. }
  305. }
  306. void FEM4EllipticPDE::SolveOLD(const std::string& fname) {
  307. Real r0[3];
  308. Real r1[3];
  309. /////////////////////////////////////////////////////////////////////////
  310. // Surface filter, to determine if points are on boundary, and need
  311. // boundary conditions applied
  312. vtkDataSetSurfaceFilter* Surface = vtkDataSetSurfaceFilter::New();
  313. Surface->SetInputData(vtkGrid);
  314. Surface->PassThroughPointIdsOn( );
  315. Surface->Update();
  316. vtkIdTypeArray* BdryIds = static_cast<vtkIdTypeArray*>
  317. (Surface->GetOutput()->GetPointData()->GetScalars("vtkOriginalPointIds"));
  318. // Expensive search for whether or not point is on boundary. O(n) cost.
  319. VectorXi bndryCnt = VectorXi::Zero(vtkGrid->GetNumberOfPoints());
  320. for (int isp=0; isp < Surface->GetOutput()->GetNumberOfPoints(); ++isp) {
  321. //double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  322. // x \in -14.5 to 14.5
  323. // y \in 0 to 30
  324. bndryCnt(BdryIds->GetTuple1(isp)) += 1;
  325. }
  326. /////////////////////////////////////////////////////////////////////////
  327. // Build stiffness matrix (A)
  328. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  329. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  330. Eigen::SparseMatrix<Real> A(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  331. std::vector< Eigen::Triplet<Real> > coeffs;
  332. // Here we iterate over all of the cells
  333. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  334. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  335. // TODO, in production code we might not want to do this check here
  336. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  337. std::cout << "DOOM FEM4EllipticPDE encountered non-tetrahedral mesh\n";
  338. std::cout << "Number of points in cell " << vtkGrid->GetCell(ic)->GetNumberOfPoints() << std::endl ;
  339. exit(1);
  340. }
  341. // construct coordinate matrix C
  342. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  343. for (int ip=0; ip<4; ++ip) {
  344. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  345. C(ip, 0) = 1;
  346. C(ip, 1) = pts[0] ;
  347. C(ip, 2) = pts[1] ;
  348. C(ip, 3) = pts[2] ;
  349. }
  350. Eigen::Matrix<Real, 4, 4> Phi = C.inverse(); // nabla \phi
  351. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  352. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  353. int ID[4];
  354. ID[0] = Ids->GetId(0);
  355. ID[1] = Ids->GetId(1);
  356. ID[2] = Ids->GetId(2);
  357. ID[3] = Ids->GetId(3);
  358. Real sum(0);
  359. Real sigma_bar = vtkGrid->GetCellData()->GetScalars()->GetTuple1(ic);
  360. for (int ip=0; ip<4; ++ip) {
  361. for (int ip2=0; ip2<4; ++ip2) {
  362. if (ip2 == ip) {
  363. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  364. // solve for the boundaries? Is one better? This seems to work, which is nice.
  365. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ip] ); // + sum;
  366. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ip])[0];
  367. //std::cout << "bb " << bb << std::endl;
  368. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  369. coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], bdry + Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  370. } else {
  371. coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  372. }
  373. // Stiffness matrix no longer contains boundary conditions...
  374. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ip], ID[ip2], Phi.col(ip).tail<3>().dot(Phi.col(ip2).tail<3>() ) * V * sigma_bar ) );
  375. }
  376. }
  377. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  378. }
  379. A.setFromTriplets(coeffs.begin(), coeffs.end());
  380. //A.makeCompressed();
  381. ////////////////////////////////////////////////////////////
  382. // Load vector g, solution vector u
  383. std::cout << "\nBuilding load vector (g)" << std::endl;
  384. VectorXr g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  385. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  386. // If the G function has been evaluated at each *node*
  387. // --> but still needs to be integrated at the surfaces
  388. // Aha, requires that there is in fact a pointdata memeber // BUG TODO BUG!!!
  389. std::cout << "Point Data ptr " << vtkGrid->GetPointData() << std::endl;
  390. //if ( vtkGrid->GetPointData() != NULL && std::string( vtkGrid->GetPointData()->GetScalars()->GetName() ).compare( std::string("G") ) == 0 ) {
  391. bool pe(false);
  392. bool ne(false);
  393. if ( true ) {
  394. std::cout << "\nUsing G from file" << std::endl;
  395. /* 3D Phi */
  396. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  397. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  398. for (int ip=0; ip<4; ++ip) {
  399. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ip); //[ipc] ;
  400. C(ip, 0) = 1;
  401. C(ip, 1) = pts[0];
  402. C(ip, 2) = pts[1];
  403. C(ip, 3) = pts[2];
  404. }
  405. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  406. //Real V = C.determinant(); // volume of tetrahedra
  407. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  408. int ID[4];
  409. ID[0] = Ids->GetId(0);
  410. ID[1] = Ids->GetId(1);
  411. ID[2] = Ids->GetId(2);
  412. ID[3] = Ids->GetId(3);
  413. /* bad news bears for magnet */
  414. double* pt = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(0);
  415. Real sum(0);
  416. /*
  417. if (!pe) {
  418. if (std::abs(pt[0]) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  419. sum = 1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0];
  420. pe = true;
  421. }
  422. }*/
  423. if (ID[0] == 26) {
  424. sum = 10;
  425. }
  426. if (ID[0] == 30) {
  427. sum = -10;
  428. }
  429. /*
  430. if (!ne) {
  431. if (std::abs(pt[0]+1.) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  432. sum = -1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0];
  433. std::cout << "Negative Electroce\n";
  434. ne = true;
  435. }
  436. }
  437. */
  438. //for (int ip=0; ip<4; ++ip) {
  439. //g(ID[ip]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0] ) ;
  440. //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ip])[0]) > 1e-3 )
  441. //}
  442. // TODO check Load Vector...
  443. g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  444. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  445. }
  446. /*
  447. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  448. vtkGrid->GetPoint(ic, r0);
  449. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  450. double g0 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ;
  451. //std::cout << "num conn " << connectedVertices->GetNumberOfIds() << std::endl;
  452. if ( std::abs(g0) > 1e-3 ) {
  453. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  454. int ip = connectedVertices->GetId(i);
  455. vtkGrid->GetPoint(ip, r1);
  456. double g1 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ip)[0] ;
  457. //g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  458. if ( std::abs(g1) > 1e-3 ) {
  459. g(ic) += CompositeSimpsons2(g1, g0, r1, r0, 1000);
  460. }
  461. //g(ic) += CompositeSimpsons2(g0, r1, r0, 8);
  462. //if ( std::abs(g1) > 1e-3 ) {
  463. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 8);
  464. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 100); // / (2*dist(r0,r1)) ;
  465. // g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  466. //g(ic) += CompositeSimpsons2(g0, r0, r1, 8);
  467. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  468. //} //else {
  469. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  470. //}
  471. }
  472. }
  473. //g(ic) = 2.* vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ; // Why 2?
  474. //std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  475. }
  476. */
  477. } else if (vtkG) { // VTK implicit function, proceed with care
  478. std::cout << "\nUsing implicit file from file" << std::endl;
  479. // OpenMP right here
  480. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  481. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  482. //vtkGrid->GetPoint(ic, r0);
  483. //g(ic) += vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  484. // std::cout << vtkG->FunctionValue(r0[0], r0[1], r0[2]) << std::endl;
  485. //g(ic) += vtkGrid->GetPointData()->GetScalars()->GetTuple1(ic);// FunctionValue(r0[0], r0[1], r0[2]) ;
  486. /*
  487. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  488. int ip = connectedVertices->GetId(i);
  489. vtkGrid->GetPoint(ip, r1);
  490. g(ic) += CompositeSimpsons2(vtkG, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  491. }
  492. */
  493. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  494. }
  495. } else if (gFcn3) {
  496. std::cout << "\nUsing gFcn3" << std::endl;
  497. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  498. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  499. vtkGrid->GetPoint(ic, r0);
  500. // TODO, test OMP sum reduction here. Is vtkGrid->GetPoint thread safe?
  501. //Real sum(0.);
  502. //#ifdef LEMMAUSEOMP
  503. //#pragma omp parallel for reduction(+:sum)
  504. //#endif
  505. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  506. int ip = connectedVertices->GetId(i);
  507. vtkGrid->GetPoint(ip, r1);
  508. g(ic) += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  509. //sum += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  510. }
  511. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  512. //g(ic) = sum;
  513. }
  514. } else {
  515. std::cout << "No source specified\n";
  516. exit(EXIT_FAILURE);
  517. }
  518. // std::cout << g << std::endl;
  519. //g(85) = 1;
  520. std::cout << "\nSolving" << std::endl;
  521. ////////////////////////////////////////////////////////////
  522. // Solving:
  523. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  524. //VectorXr u = chol.solve(g);
  525. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  526. Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  527. //Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  528. cg.setMaxIterations(3000);
  529. //cg.compute(A);
  530. //std::cout << "Computed " << std::endl;
  531. VectorXr u = cg.solve(g);
  532. std::cout << "#iterations: " << cg.iterations() << std::endl;
  533. std::cout << "estimated error: " << cg.error() << std::endl;
  534. vtkDoubleArray *gArray = vtkDoubleArray::New();
  535. vtkDoubleArray *uArray = vtkDoubleArray::New();
  536. uArray->SetNumberOfComponents(1);
  537. gArray->SetNumberOfComponents(1);
  538. for (int iu = 0; iu<u.size(); ++iu) {
  539. uArray->InsertTuple1(iu, u[iu]);
  540. gArray->InsertTuple1(iu, g[iu]);
  541. }
  542. uArray->SetName("u");
  543. gArray->SetName("g");
  544. vtkGrid->GetPointData()->AddArray(uArray);
  545. vtkGrid->GetPointData()->AddArray(gArray);
  546. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  547. Writer->SetInputData(vtkGrid);
  548. Writer->SetFileName(fname.c_str());
  549. Writer->Write();
  550. Writer->Delete();
  551. Surface->Delete();
  552. gArray->Delete();
  553. uArray->Delete();
  554. }
  555. // Uses simpon's rule to perform a definite integral of a
  556. // function (passed as a pointer). The integration occurs from
  557. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  558. Real FEM4EllipticPDE::CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int n) {
  559. Vector3r R0(r0[0], r0[1], r0[2]);
  560. Vector3r R1(r1[0], r1[1], r1[2]);
  561. // force n to be even
  562. assert(n > 0);
  563. //n += (n % 2);
  564. Real h = dist(r0, r1) / (Real)(n) ;
  565. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  566. Vector3r dr = (R1 - R0).array() / Real(n);
  567. Vector3r rx;
  568. rx.array() = R0.array() + dr.array();
  569. for (int i=1; i<n; i+=2) {
  570. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]);
  571. rx += 2.*dr;
  572. }
  573. rx.array() = R0.array() + 2*dr.array();
  574. for (int i=2; i<n-1; i+=2) {
  575. S += 2.*f->FunctionValue(rx[0], rx[1], rx[2]);
  576. rx += 2.*dr;
  577. }
  578. return h * S / 3.;
  579. }
  580. // Uses simpon's rule to perform a definite integral of a
  581. // function (passed as a pointer). The integration occurs from
  582. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  583. // This is just here as a convenience
  584. Real FEM4EllipticPDE::CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int n) {
  585. return dist(r0,r1)*f;
  586. /*
  587. Vector3r R0(r0[0], r0[1], r0[2]);
  588. Vector3r R1(r1[0], r1[1], r1[2]);
  589. // force n to be even
  590. assert(n > 0);
  591. //n += (n % 2);
  592. Real h = dist(r0, r1) / (Real)(n) ;
  593. Real S = f + f;
  594. Vector3r dr = (R1 - R0).array() / Real(n);
  595. //Vector3r rx;
  596. //rx.array() = R0.array() + dr.array();
  597. for (int i=1; i<n; i+=2) {
  598. S += 4. * f;
  599. //rx += 2.*dr;
  600. }
  601. //rx.array() = R0.array() + 2*dr.array();
  602. for (int i=2; i<n-1; i+=2) {
  603. S += 2. * f;
  604. //rx += 2.*dr;
  605. }
  606. return h * S / 3.;
  607. */
  608. }
  609. /*
  610. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  611. * test function owned by the FEM implimentaion.
  612. */
  613. Real FEM4EllipticPDE::CompositeSimpsons2(vtkImplicitFunction* f,
  614. Real r0[3], Real r1[3], int n) {
  615. Vector3r R0(r0[0], r0[1], r0[2]);
  616. Vector3r R1(r1[0], r1[1], r1[2]);
  617. // force n to be even
  618. assert(n > 0);
  619. //n += (n % 2);
  620. Real h = dist(r0, r1) / (Real)(n) ;
  621. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  622. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  623. //Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  624. Vector3r dr = (R1 - R0).array() / Real(n);
  625. Vector3r rx;
  626. rx.array() = R0.array() + dr.array();
  627. for (int i=1; i<n; i+=2) {
  628. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  629. rx += 2.*dr;
  630. }
  631. rx.array() = R0.array() + 2*dr.array();
  632. for (int i=2; i<n-1; i+=2) {
  633. S += 2. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  634. rx += 2.*dr;
  635. }
  636. return h * S / 3.;
  637. }
  638. /*
  639. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  640. * test function owned by the FEM implimentaion.
  641. */
  642. Real FEM4EllipticPDE::CompositeSimpsons2( Real (*f)(const Real&, const Real&, const Real&),
  643. Real r0[3], Real r1[3], int n) {
  644. Vector3r R0(r0[0], r0[1], r0[2]);
  645. Vector3r R1(r1[0], r1[1], r1[2]);
  646. // force n to be even
  647. assert(n > 0);
  648. //n += (n % 2);
  649. Real h = dist(r0, r1) / (Real)(n) ;
  650. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  651. //Real S = f(r0[0], r0[1], r0[2])*Hat(R0, R0, R1) + f(r1[0], r1[1], r1[2])*Hat(R1, R0, R1);
  652. Real S = f(r0[0], r0[1], r0[2]) + f(r1[0], r1[1], r1[2]);
  653. Vector3r dr = (R1 - R0).array() / Real(n);
  654. Vector3r rx;
  655. rx.array() = R0.array() + dr.array();
  656. for (int i=1; i<n; i+=2) {
  657. S += 4. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  658. rx += 2.*dr;
  659. }
  660. rx.array() = R0.array() + 2*dr.array();
  661. for (int i=2; i<n-1; i+=2) {
  662. S += 2. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  663. rx += 2.*dr;
  664. }
  665. return h * S / 3.;
  666. }
  667. /*
  668. * Performs numerical integration of two functions, one is constant valued f, the other is the FEM
  669. * test function owned by the FEM implimentaion.
  670. */
  671. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int n) {
  672. Vector3r R0(r0[0], r0[1], r0[2]);
  673. Vector3r R1(r1[0], r1[1], r1[2]);
  674. // force n to be even
  675. assert(n > 0);
  676. //n += (n % 2);
  677. Real h = dist(r0, r1) / (Real)(n) ;
  678. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  679. Real S = 2*f; //*Hat(R0, R0, R1) + f*Hat(R1, R0, R1);
  680. Vector3r dr = (R1 - R0).array() / Real(n);
  681. Vector3r rx;
  682. rx.array() = R0.array() + dr.array();
  683. for (int i=1; i<n; i+=2) {
  684. S += 4. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  685. rx += 2.*dr;
  686. }
  687. rx.array() = R0.array() + 2*dr.array();
  688. for (int i=2; i<n-1; i+=2) {
  689. S += 2. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  690. rx += 2.*dr;
  691. }
  692. return h * S / 3.;
  693. }
  694. /*
  695. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  696. * test function owned by the FEM implimentaion.
  697. */
  698. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  699. Vector3r R0(r0[0], r0[1], r0[2]);
  700. Vector3r R1(r1[0], r1[1], r1[2]);
  701. // force n to be even
  702. assert(n > 0);
  703. //n += (n % 2);
  704. Real h = dist(r0, r1) / (Real)(n) ;
  705. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  706. // NO! We are looking at 1/2 hat often!!! So S = f1!
  707. Real S = f1; //f0*Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  708. Vector3r dr = (R1 - R0).array() / Real(n);
  709. // linear interpolate function
  710. //Vector3r rx;
  711. //rx.array() = R0.array() + dr.array();
  712. for (int i=1; i<n; i+=2) {
  713. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  714. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0) ;
  715. }
  716. //rx.array() = R0.array() + 2*dr.array();
  717. for (int i=2; i<n-1; i+=2) {
  718. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  719. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0);
  720. }
  721. return h * S / 3.;
  722. }
  723. /*
  724. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  725. * test function owned by the FEM implimentaion.
  726. */
  727. Real FEM4EllipticPDE::CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  728. Vector3r R0(r0[0], r0[1], r0[2]);
  729. Vector3r R1(r1[0], r1[1], r1[2]);
  730. // force n to be even
  731. assert(n > 0);
  732. //n += (n % 2);
  733. Real h = dist(r0, r1) / (Real)(n) ;
  734. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  735. // NO! We are looking at 1/2 hat often!!! So S = f1!
  736. Real S = f0+f1; //Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  737. Vector3r dr = (R1 - R0).array() / Real(n);
  738. // linear interpolate function
  739. //Vector3r rx;
  740. //rx.array() = R0.array() + dr.array();
  741. for (int i=1; i<n; i+=2) {
  742. double fx = 1;// f0 + (f1 - f0) * ((i*h)/(h*n));
  743. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1) * Hat(R1.array() + i*dr.array(), R1, R0) ;
  744. }
  745. //rx.array() = R0.array() + 2*dr.array();
  746. for (int i=2; i<n-1; i+=2) {
  747. double fx = 1; // f0 + (f1 - f0) * ((i*h)/(h*n));
  748. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1)* Hat(R1.array() + i*dr.array(), R1, R0);
  749. }
  750. return h * S / 3.;
  751. }
  752. //--------------------------------------------------------------------------------------
  753. // Class: FEM4EllipticPDE
  754. // Method: Hat
  755. //--------------------------------------------------------------------------------------
  756. Real FEM4EllipticPDE::Hat ( const Vector3r& r, const Vector3r& r0, const Vector3r& r1 ) {
  757. //return (r-r0).norm() / (r1-r0).norm() ;
  758. return dist(r, r0) / dist(r1, r0) ;
  759. } // ----- end of method FEM4EllipticPDE::Hat -----
  760. } // ----- end of Lemma name -----