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

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