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

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