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 38KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954
  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 ii=0; ii<4; ++ii) {
  131. // double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  132. // C(ii, 0) = 1;
  133. // C(ii, 1) = pts[0];
  134. // C(ii, 2) = pts[1];
  135. // C(ii, 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 <<" rows\tcols\n";
  177. std::cout << "A: " << A.rows() << "\t" << A.cols() << std::endl;
  178. std::cout << "g: " << g.rows() << "\t" << g.cols() << std::endl;
  179. VectorXr u = cg.solve(g);
  180. std::cout << "#iterations: " << cg.iterations() << std::endl;
  181. std::cout << "estimated error: " << cg.error() << std::endl;
  182. vtkDoubleArray *gArray = vtkDoubleArray::New();
  183. vtkDoubleArray *uArray = vtkDoubleArray::New();
  184. uArray->SetNumberOfComponents(1);
  185. gArray->SetNumberOfComponents(1);
  186. for (int iu = 0; iu<u.size(); ++iu) {
  187. uArray->InsertTuple1(iu, u[iu]);
  188. gArray->InsertTuple1(iu, g[iu]);
  189. }
  190. uArray->SetName("u");
  191. gArray->SetName("g");
  192. vtkGrid->GetPointData()->AddArray(uArray);
  193. vtkGrid->GetPointData()->AddArray(gArray);
  194. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  195. Writer->SetInputData(vtkGrid);
  196. Writer->SetFileName(resfile.c_str());
  197. Writer->Write();
  198. Writer->Delete();
  199. gArray->Delete();
  200. uArray->Delete();
  201. }
  202. //--------------------------------------------------------------------------------------
  203. // Class: FEM4EllipticPDE
  204. // Method: ConstructAMatrix
  205. //--------------------------------------------------------------------------------------
  206. void FEM4EllipticPDE::ConstructAMatrix ( ) {
  207. /////////////////////////////////////////////////////////////////////////
  208. // Build stiffness matrix (A)
  209. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  210. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  211. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  212. //Eigen::SparseMatrix<Real>
  213. A.resize(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  214. std::vector< Eigen::Triplet<Real> > coeffs;
  215. if ( !vtkGrid->GetPointData()->GetScalars("vtkValidPointMask") ) {
  216. throw std::runtime_error("No vtkValidPointMask");
  217. }
  218. if ( !vtkGrid->GetCellData()->GetScalars("G") && !vtkGrid->GetPointData()->GetScalars("G") ) {
  219. throw std::runtime_error("No Cell or Point Data G");
  220. }
  221. bool GCell = false;
  222. if ( vtkGrid->GetCellData()->GetScalars("G") ) {
  223. GCell = true;
  224. }
  225. // Here we iterate over all of the cells
  226. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  227. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  228. // TODO, in production code we might not want to do this check here
  229. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  230. throw std::runtime_error("Non-tetrahedral mesh encountered!");
  231. }
  232. // construct coordinate matrix C
  233. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  234. for (int ii=0; ii<4; ++ii) {
  235. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  236. C(ii, 0) = 1;
  237. C(ii, 1) = pts[0] ;
  238. C(ii, 2) = pts[1] ;
  239. C(ii, 3) = pts[2] ;
  240. }
  241. Eigen::Matrix<Real, 4, 4> Phi = C.inverse(); // nabla \phi
  242. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  243. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  244. int ID[4];
  245. ID[0] = Ids->GetId(0);
  246. ID[1] = Ids->GetId(1);
  247. ID[2] = Ids->GetId(2);
  248. ID[3] = Ids->GetId(3);
  249. Real sigma_bar(0);
  250. if (GCell) {
  251. sigma_bar = vtkGrid->GetCellData()->GetScalars("G")->GetTuple1(ic);
  252. } else {
  253. sigma_bar = vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0]);
  254. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1]);
  255. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2]);
  256. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3]);
  257. sigma_bar /= 4.;
  258. }
  259. sigma_bar = 1.;
  260. for (int ii=0; ii<4; ++ii) {
  261. for (int jj=0; jj<4; ++jj) {
  262. if (jj == ii) {
  263. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  264. // solve for the boundaries? Is one better? This seems to work, which is nice.
  265. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ii] ); // + sum;
  266. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ii])[0];
  267. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  268. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], bdry + Phi.col(ii).tail<3>().dot(Phi.col(jj).tail<3>() ) * V * sigma_bar ) );
  269. } else {
  270. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], Phi.col(ii).tail<3>().dot(Phi.col(jj).tail<3>() ) * V * sigma_bar ) );
  271. }
  272. // Stiffness matrix no longer contains boundary conditions...
  273. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], Phi.col(ii).tail<3>().dot(Phi.col(jj).tail<3>() ) * V * sigma_bar ) );
  274. }
  275. }
  276. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  277. }
  278. A.setFromTriplets(coeffs.begin(), coeffs.end());
  279. //std::cout << "A\n" << A << std::endl;
  280. }
  281. void FEM4EllipticPDE::SetupPotential() {
  282. ////////////////////////////////////////////////////////////
  283. // Load vector g
  284. std::cout << "\nBuilding load vector (g)" << std::endl;
  285. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  286. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " points" << std::endl;
  287. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  288. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  289. for (int ii=0; ii<4; ++ii) {
  290. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  291. C(ii, 0) = 1;
  292. C(ii, 1) = pts[0];
  293. C(ii, 2) = pts[1];
  294. C(ii, 3) = pts[2];
  295. }
  296. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  297. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  298. int ID[4];
  299. ID[0] = Ids->GetId(0);
  300. ID[1] = Ids->GetId(1);
  301. ID[2] = Ids->GetId(2);
  302. ID[3] = Ids->GetId(3);
  303. Real avg(0);
  304. Real GG[4];
  305. for (int ii=0; ii<4; ++ii) {
  306. GG[ii] = vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0];
  307. //avg /= 4.;
  308. }
  309. if ( std::abs( (GG[0]+GG[1]+GG[2]+GG[3])/4. - GG[0]) < 1e-5) {
  310. avg = GG[0];
  311. }
  312. for (int ii=0; ii<4; ++ii) {
  313. // avg += vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  314. // //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0]) > 1e-3 )
  315. //}
  316. //TODO this seems wrong!
  317. //avg /= 4.;
  318. //g(ID[ii]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0] ) ;
  319. g(ID[ii]) += (V/4.) * avg;
  320. }
  321. }
  322. }
  323. void FEM4EllipticPDE::SolveOLD(const std::string& fname) {
  324. Real r0[3];
  325. Real r1[3];
  326. /////////////////////////////////////////////////////////////////////////
  327. // Surface filter, to determine if points are on boundary, and need
  328. // boundary conditions applied
  329. vtkDataSetSurfaceFilter* Surface = vtkDataSetSurfaceFilter::New();
  330. Surface->SetInputData(vtkGrid);
  331. Surface->PassThroughPointIdsOn( );
  332. Surface->Update();
  333. vtkIdTypeArray* BdryIds = static_cast<vtkIdTypeArray*>
  334. (Surface->GetOutput()->GetPointData()->GetScalars("vtkOriginalPointIds"));
  335. // Expensive search for whether or not point is on boundary. O(n) cost.
  336. VectorXi bndryCnt = VectorXi::Zero(vtkGrid->GetNumberOfPoints());
  337. for (int isp=0; isp < Surface->GetOutput()->GetNumberOfPoints(); ++isp) {
  338. //double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  339. // x \in -14.5 to 14.5
  340. // y \in 0 to 30
  341. bndryCnt(BdryIds->GetTuple1(isp)) += 1;
  342. }
  343. /////////////////////////////////////////////////////////////////////////
  344. // Build stiffness matrix (A)
  345. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  346. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  347. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  348. Eigen::SparseMatrix<Real> A(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  349. std::vector< Eigen::Triplet<Real> > coeffs;
  350. // Here we iterate over all of the cells
  351. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  352. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  353. // TODO, in production code we might not want to do this check here
  354. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  355. std::cout << "DOOM FEM4EllipticPDE encountered non-tetrahedral mesh\n";
  356. std::cout << "Number of points in cell " << vtkGrid->GetCell(ic)->GetNumberOfPoints() << std::endl ;
  357. exit(1);
  358. }
  359. // construct coordinate matrix C
  360. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  361. for (int ii=0; ii<4; ++ii) {
  362. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  363. C(ii, 0) = 1;
  364. C(ii, 1) = pts[0] ;
  365. C(ii, 2) = pts[1] ;
  366. C(ii, 3) = pts[2] ;
  367. }
  368. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  369. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  370. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  371. int ID[4];
  372. ID[0] = Ids->GetId(0);
  373. ID[1] = Ids->GetId(1);
  374. ID[2] = Ids->GetId(2);
  375. ID[3] = Ids->GetId(3);
  376. Real sum(0);
  377. Real sigma_bar = vtkGrid->GetCellData()->GetScalars()->GetTuple1(ic);
  378. for (int ii=0; ii<4; ++ii) {
  379. for (int jj=0; jj<4; ++jj) {
  380. if (jj == ii) {
  381. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  382. // solve for the boundaries? Is one better? This seems to work, which is nice.
  383. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ii] ); // + sum;
  384. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ii])[0];
  385. //std::cout << "bb " << bb << std::endl;
  386. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  387. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], bdry + GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  388. } else {
  389. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  390. }
  391. // Stiffness matrix no longer contains boundary conditions...
  392. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  393. }
  394. }
  395. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  396. }
  397. A.setFromTriplets(coeffs.begin(), coeffs.end());
  398. //A.makeCompressed();
  399. ////////////////////////////////////////////////////////////
  400. // Load vector g, solution vector u
  401. std::cout << "\nBuilding load vector (g)" << std::endl;
  402. VectorXr g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  403. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  404. // If the G function has been evaluated at each *node*
  405. // --> but still needs to be integrated at the surfaces
  406. // Aha, requires that there is in fact a pointdata memeber // BUG TODO BUG!!!
  407. std::cout << "Point Data ptr " << vtkGrid->GetPointData() << std::endl;
  408. //if ( vtkGrid->GetPointData() != NULL && std::string( vtkGrid->GetPointData()->GetScalars()->GetName() ).compare( std::string("G") ) == 0 ) {
  409. bool pe(false);
  410. bool ne(false);
  411. if ( true ) {
  412. std::cout << "\nUsing G from file" << std::endl;
  413. /* 3D Phi */
  414. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  415. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  416. for (int ii=0; ii<4; ++ii) {
  417. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  418. C(ii, 0) = 1;
  419. C(ii, 1) = pts[0];
  420. C(ii, 2) = pts[1];
  421. C(ii, 3) = pts[2];
  422. }
  423. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  424. //Real V = C.determinant(); // volume of tetrahedra
  425. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  426. int ID[4];
  427. ID[0] = Ids->GetId(0);
  428. ID[1] = Ids->GetId(1);
  429. ID[2] = Ids->GetId(2);
  430. ID[3] = Ids->GetId(3);
  431. /* bad news bears for magnet */
  432. double* pt = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(0);
  433. Real sum(0);
  434. /*
  435. if (!pe) {
  436. if (std::abs(pt[0]) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  437. sum = 1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  438. pe = true;
  439. }
  440. }*/
  441. if (ID[0] == 26) {
  442. sum = 10;
  443. }
  444. if (ID[0] == 30) {
  445. sum = -10;
  446. }
  447. /*
  448. if (!ne) {
  449. if (std::abs(pt[0]+1.) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  450. sum = -1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  451. std::cout << "Negative Electroce\n";
  452. ne = true;
  453. }
  454. }
  455. */
  456. //for (int ii=0; ii<4; ++ii) {
  457. //g(ID[ii]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0] ) ;
  458. //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0]) > 1e-3 )
  459. //}
  460. // TODO check Load Vector...
  461. g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  462. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  463. }
  464. /*
  465. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  466. vtkGrid->GetPoint(ic, r0);
  467. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  468. double g0 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ;
  469. //std::cout << "num conn " << connectedVertices->GetNumberOfIds() << std::endl;
  470. if ( std::abs(g0) > 1e-3 ) {
  471. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  472. int ii = connectedVertices->GetId(i);
  473. vtkGrid->GetPoint(ii, r1);
  474. double g1 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ii)[0] ;
  475. //g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  476. if ( std::abs(g1) > 1e-3 ) {
  477. g(ic) += CompositeSimpsons2(g1, g0, r1, r0, 1000);
  478. }
  479. //g(ic) += CompositeSimpsons2(g0, r1, r0, 8);
  480. //if ( std::abs(g1) > 1e-3 ) {
  481. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 8);
  482. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 100); // / (2*dist(r0,r1)) ;
  483. // g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  484. //g(ic) += CompositeSimpsons2(g0, r0, r1, 8);
  485. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  486. //} //else {
  487. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  488. //}
  489. }
  490. }
  491. //g(ic) = 2.* vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ; // Why 2?
  492. //std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  493. }
  494. */
  495. } else if (vtkG) { // VTK implicit function, proceed with care
  496. std::cout << "\nUsing implicit file from file" << std::endl;
  497. // OpenMP right here
  498. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  499. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  500. //vtkGrid->GetPoint(ic, r0);
  501. //g(ic) += vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  502. // std::cout << vtkG->FunctionValue(r0[0], r0[1], r0[2]) << std::endl;
  503. //g(ic) += vtkGrid->GetPointData()->GetScalars()->GetTuple1(ic);// FunctionValue(r0[0], r0[1], r0[2]) ;
  504. /*
  505. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  506. int ii = connectedVertices->GetId(i);
  507. vtkGrid->GetPoint(ii, r1);
  508. g(ic) += CompositeSimpsons2(vtkG, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  509. }
  510. */
  511. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  512. }
  513. } else if (gFcn3) {
  514. std::cout << "\nUsing gFcn3" << std::endl;
  515. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  516. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  517. vtkGrid->GetPoint(ic, r0);
  518. // TODO, test OMP sum reduction here. Is vtkGrid->GetPoint thread safe?
  519. //Real sum(0.);
  520. //#ifdef LEMMAUSEOMP
  521. //#pragma omp parallel for reduction(+:sum)
  522. //#endif
  523. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  524. int ii = connectedVertices->GetId(i);
  525. vtkGrid->GetPoint(ii, r1);
  526. g(ic) += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  527. //sum += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  528. }
  529. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  530. //g(ic) = sum;
  531. }
  532. } else {
  533. std::cout << "No source specified\n";
  534. exit(EXIT_FAILURE);
  535. }
  536. // std::cout << g << std::endl;
  537. //g(85) = 1;
  538. std::cout << "\nSolving" << std::endl;
  539. ////////////////////////////////////////////////////////////
  540. // Solving:
  541. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  542. //VectorXr u = chol.solve(g);
  543. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  544. Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  545. //Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  546. cg.setMaxIterations(3000);
  547. //cg.compute(A);
  548. //std::cout << "Computed " << std::endl;
  549. VectorXr u = cg.solve(g);
  550. std::cout << "#iterations: " << cg.iterations() << std::endl;
  551. std::cout << "estimated error: " << cg.error() << std::endl;
  552. vtkDoubleArray *gArray = vtkDoubleArray::New();
  553. vtkDoubleArray *uArray = vtkDoubleArray::New();
  554. uArray->SetNumberOfComponents(1);
  555. gArray->SetNumberOfComponents(1);
  556. for (int iu = 0; iu<u.size(); ++iu) {
  557. uArray->InsertTuple1(iu, u[iu]);
  558. gArray->InsertTuple1(iu, g[iu]);
  559. }
  560. uArray->SetName("u");
  561. gArray->SetName("g");
  562. vtkGrid->GetPointData()->AddArray(uArray);
  563. vtkGrid->GetPointData()->AddArray(gArray);
  564. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  565. Writer->SetInputData(vtkGrid);
  566. Writer->SetFileName(fname.c_str());
  567. Writer->Write();
  568. Writer->Delete();
  569. Surface->Delete();
  570. gArray->Delete();
  571. uArray->Delete();
  572. }
  573. // Uses simpon's rule to perform a definite integral of a
  574. // function (passed as a pointer). The integration occurs from
  575. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  576. Real FEM4EllipticPDE::CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int n) {
  577. Vector3r R0(r0[0], r0[1], r0[2]);
  578. Vector3r R1(r1[0], r1[1], r1[2]);
  579. // force n to be even
  580. assert(n > 0);
  581. //n += (n % 2);
  582. Real h = dist(r0, r1) / (Real)(n) ;
  583. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  584. Vector3r dr = (R1 - R0).array() / Real(n);
  585. Vector3r rx;
  586. rx.array() = R0.array() + dr.array();
  587. for (int i=1; i<n; i+=2) {
  588. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]);
  589. rx += 2.*dr;
  590. }
  591. rx.array() = R0.array() + 2*dr.array();
  592. for (int i=2; i<n-1; i+=2) {
  593. S += 2.*f->FunctionValue(rx[0], rx[1], rx[2]);
  594. rx += 2.*dr;
  595. }
  596. return h * S / 3.;
  597. }
  598. // Uses simpon's rule to perform a definite integral of a
  599. // function (passed as a pointer). The integration occurs from
  600. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  601. // This is just here as a convenience
  602. Real FEM4EllipticPDE::CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int n) {
  603. return dist(r0,r1)*f;
  604. /*
  605. Vector3r R0(r0[0], r0[1], r0[2]);
  606. Vector3r R1(r1[0], r1[1], r1[2]);
  607. // force n to be even
  608. assert(n > 0);
  609. //n += (n % 2);
  610. Real h = dist(r0, r1) / (Real)(n) ;
  611. Real S = f + f;
  612. Vector3r dr = (R1 - R0).array() / Real(n);
  613. //Vector3r rx;
  614. //rx.array() = R0.array() + dr.array();
  615. for (int i=1; i<n; i+=2) {
  616. S += 4. * f;
  617. //rx += 2.*dr;
  618. }
  619. //rx.array() = R0.array() + 2*dr.array();
  620. for (int i=2; i<n-1; i+=2) {
  621. S += 2. * f;
  622. //rx += 2.*dr;
  623. }
  624. return h * S / 3.;
  625. */
  626. }
  627. /*
  628. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  629. * test function owned by the FEM implimentaion.
  630. */
  631. Real FEM4EllipticPDE::CompositeSimpsons2(vtkImplicitFunction* f,
  632. Real r0[3], Real r1[3], int n) {
  633. Vector3r R0(r0[0], r0[1], r0[2]);
  634. Vector3r R1(r1[0], r1[1], r1[2]);
  635. // force n to be even
  636. assert(n > 0);
  637. //n += (n % 2);
  638. Real h = dist(r0, r1) / (Real)(n) ;
  639. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  640. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  641. //Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  642. Vector3r dr = (R1 - R0).array() / Real(n);
  643. Vector3r rx;
  644. rx.array() = R0.array() + dr.array();
  645. for (int i=1; i<n; i+=2) {
  646. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  647. rx += 2.*dr;
  648. }
  649. rx.array() = R0.array() + 2*dr.array();
  650. for (int i=2; i<n-1; i+=2) {
  651. S += 2. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  652. rx += 2.*dr;
  653. }
  654. return h * S / 3.;
  655. }
  656. /*
  657. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  658. * test function owned by the FEM implimentaion.
  659. */
  660. Real FEM4EllipticPDE::CompositeSimpsons2( Real (*f)(const Real&, const Real&, const Real&),
  661. Real r0[3], Real r1[3], int n) {
  662. Vector3r R0(r0[0], r0[1], r0[2]);
  663. Vector3r R1(r1[0], r1[1], r1[2]);
  664. // force n to be even
  665. assert(n > 0);
  666. //n += (n % 2);
  667. Real h = dist(r0, r1) / (Real)(n) ;
  668. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  669. //Real S = f(r0[0], r0[1], r0[2])*Hat(R0, R0, R1) + f(r1[0], r1[1], r1[2])*Hat(R1, R0, R1);
  670. Real S = f(r0[0], r0[1], r0[2]) + f(r1[0], r1[1], r1[2]);
  671. Vector3r dr = (R1 - R0).array() / Real(n);
  672. Vector3r rx;
  673. rx.array() = R0.array() + dr.array();
  674. for (int i=1; i<n; i+=2) {
  675. S += 4. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  676. rx += 2.*dr;
  677. }
  678. rx.array() = R0.array() + 2*dr.array();
  679. for (int i=2; i<n-1; i+=2) {
  680. S += 2. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  681. rx += 2.*dr;
  682. }
  683. return h * S / 3.;
  684. }
  685. /*
  686. * Performs numerical integration of two functions, one is constant valued f, the other is the FEM
  687. * test function owned by the FEM implimentaion.
  688. */
  689. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int n) {
  690. Vector3r R0(r0[0], r0[1], r0[2]);
  691. Vector3r R1(r1[0], r1[1], r1[2]);
  692. // force n to be even
  693. assert(n > 0);
  694. //n += (n % 2);
  695. Real h = dist(r0, r1) / (Real)(n) ;
  696. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  697. Real S = 2*f; //*Hat(R0, R0, R1) + f*Hat(R1, R0, R1);
  698. Vector3r dr = (R1 - R0).array() / Real(n);
  699. Vector3r rx;
  700. rx.array() = R0.array() + dr.array();
  701. for (int i=1; i<n; i+=2) {
  702. S += 4. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  703. rx += 2.*dr;
  704. }
  705. rx.array() = R0.array() + 2*dr.array();
  706. for (int i=2; i<n-1; i+=2) {
  707. S += 2. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  708. rx += 2.*dr;
  709. }
  710. return h * S / 3.;
  711. }
  712. /*
  713. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  714. * test function owned by the FEM implimentaion.
  715. */
  716. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  717. Vector3r R0(r0[0], r0[1], r0[2]);
  718. Vector3r R1(r1[0], r1[1], r1[2]);
  719. // force n to be even
  720. assert(n > 0);
  721. //n += (n % 2);
  722. Real h = dist(r0, r1) / (Real)(n) ;
  723. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  724. // NO! We are looking at 1/2 hat often!!! So S = f1!
  725. Real S = f1; //f0*Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  726. Vector3r dr = (R1 - R0).array() / Real(n);
  727. // linear interpolate function
  728. //Vector3r rx;
  729. //rx.array() = R0.array() + dr.array();
  730. for (int i=1; i<n; i+=2) {
  731. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  732. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0) ;
  733. }
  734. //rx.array() = R0.array() + 2*dr.array();
  735. for (int i=2; i<n-1; i+=2) {
  736. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  737. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0);
  738. }
  739. return h * S / 3.;
  740. }
  741. /*
  742. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  743. * test function owned by the FEM implimentaion.
  744. */
  745. Real FEM4EllipticPDE::CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  746. Vector3r R0(r0[0], r0[1], r0[2]);
  747. Vector3r R1(r1[0], r1[1], r1[2]);
  748. // force n to be even
  749. assert(n > 0);
  750. //n += (n % 2);
  751. Real h = dist(r0, r1) / (Real)(n) ;
  752. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  753. // NO! We are looking at 1/2 hat often!!! So S = f1!
  754. Real S = f0+f1; //Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  755. Vector3r dr = (R1 - R0).array() / Real(n);
  756. // linear interpolate function
  757. //Vector3r rx;
  758. //rx.array() = R0.array() + dr.array();
  759. for (int i=1; i<n; i+=2) {
  760. double fx = 1;// f0 + (f1 - f0) * ((i*h)/(h*n));
  761. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1) * Hat(R1.array() + i*dr.array(), R1, R0) ;
  762. }
  763. //rx.array() = R0.array() + 2*dr.array();
  764. for (int i=2; i<n-1; i+=2) {
  765. double fx = 1; // f0 + (f1 - f0) * ((i*h)/(h*n));
  766. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1)* Hat(R1.array() + i*dr.array(), R1, R0);
  767. }
  768. return h * S / 3.;
  769. }
  770. //--------------------------------------------------------------------------------------
  771. // Class: FEM4EllipticPDE
  772. // Method: Hat
  773. //--------------------------------------------------------------------------------------
  774. Real FEM4EllipticPDE::Hat ( const Vector3r& r, const Vector3r& r0, const Vector3r& r1 ) {
  775. //return (r-r0).norm() / (r1-r0).norm() ;
  776. return dist(r, r0) / dist(r1, r0) ;
  777. } // ----- end of method FEM4EllipticPDE::Hat -----
  778. } // ----- end of Lemma name -----