Galerkin FEM for elliptic PDEs
Vous ne pouvez pas sélectionner plus de 25 sujets Les noms de sujets doivent commencer par une lettre ou un nombre, peuvent contenir des tirets ('-') et peuvent comporter jusqu'à 35 caractères.

FEM4EllipticPDE.cpp 38KB

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