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

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201
  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. #ifdef HAVE_YAMLCPP
  38. std::ostream &operator << (std::ostream &stream, const FEM4EllipticPDE &ob) {
  39. stream << ob.Serialize() << "\n---\n"; // End of doc --- as a direct stream should encapulste thingy
  40. return stream;
  41. }
  42. #else
  43. std::ostream &operator<<(std::ostream &stream,
  44. const FEM4EllipticPDE &ob) {
  45. stream << *(LemmaObject*)(&ob);
  46. return stream;
  47. }
  48. #endif
  49. // ==================== LIFECYCLE =======================
  50. FEM4EllipticPDE::FEM4EllipticPDE(const std::string&name) :
  51. LemmaObject(name), BndryH(10), BndrySigma(10),
  52. vtkSigma(nullptr), vtkG(nullptr), vtkGrid(nullptr), gFcn3(nullptr) {
  53. }
  54. #ifdef HAVE_YAMLCPP
  55. //--------------------------------------------------------------------------------------
  56. // Class: FEM4EllipticPDE
  57. // Method: FEM4EllipticPDE
  58. // Description: DeSerializing constructor (protected)
  59. //--------------------------------------------------------------------------------------
  60. FEM4EllipticPDE::FEM4EllipticPDE (const YAML::Node& node) : LemmaObject(node) {
  61. // TODO impliment
  62. throw std::runtime_error("FEM4ELLIPTICPDE YAML CONSTRUCTOR NOT IMPLIMENTED!");
  63. } // ----- end of method FEM4EllipticPDE::FEM4EllipticPDE (constructor) -----
  64. #endif
  65. FEM4EllipticPDE::~FEM4EllipticPDE() {
  66. }
  67. void FEM4EllipticPDE::Release() {
  68. delete this;
  69. }
  70. FEM4EllipticPDE* FEM4EllipticPDE::New( ) {
  71. FEM4EllipticPDE* Obj = new FEM4EllipticPDE("FEM4EllipticPDE");
  72. Obj->AttachTo(Obj);
  73. return Obj;
  74. }
  75. void FEM4EllipticPDE::Delete() {
  76. this->DetachFrom(this);
  77. }
  78. #ifdef HAVE_YAMLCPP
  79. //--------------------------------------------------------------------------------------
  80. // Class: FEM4EllipticPDE
  81. // Method: Serialize
  82. //--------------------------------------------------------------------------------------
  83. YAML::Node FEM4EllipticPDE::Serialize ( ) const {
  84. YAML::Node node = LemmaObject::Serialize();;
  85. node.SetTag( this->Name );
  86. // FILL IN CLASS SPECIFICS HERE
  87. return node;
  88. } // ----- end of method FEM4EllipticPDE::Serialize -----
  89. //--------------------------------------------------------------------------------------
  90. // Class: FEM4EllipticPDE
  91. // Method: DeSerialize
  92. //--------------------------------------------------------------------------------------
  93. FEM4EllipticPDE* FEM4EllipticPDE::DeSerialize ( const YAML::Node& node ) {
  94. FEM4EllipticPDE* Object = new FEM4EllipticPDE(node);
  95. Object->AttachTo(Object);
  96. DESERIALIZECHECK( node, Object )
  97. return Object ;
  98. } // ----- end of method FEM4EllipticPDE::DeSerialize -----
  99. #endif
  100. // ==================== OPERATIONS =======================
  101. void FEM4EllipticPDE::SetSigmaFunction(vtkImplicitFunction* sigma) {
  102. vtkSigma = sigma;
  103. }
  104. void FEM4EllipticPDE::SetBoundaryStep(const Real& h) {
  105. BndryH = h;
  106. }
  107. void FEM4EllipticPDE::SetGFunction(vtkImplicitFunction* g) {
  108. vtkG = g;
  109. }
  110. void FEM4EllipticPDE::SetGFunction( Real (*gFcn)(const Real&, const Real&, const Real&) ) {
  111. // vtkG = g;
  112. gFcn3 = gFcn;
  113. }
  114. //void FEM4EllipticPDE::SetGrid(vtkDataSet* grid) {
  115. // vtkGrid = grid;
  116. //}
  117. void FEM4EllipticPDE::SetGrid(vtkUnstructuredGrid* grid) {
  118. vtkGrid = grid;
  119. }
  120. //--------------------------------------------------------------------------------------
  121. // Class: FEM4EllipticPDE
  122. // Method: SetVTUGridFile
  123. //--------------------------------------------------------------------------------------
  124. void FEM4EllipticPDE::SetVTUGridFile ( std::string& fname ) {
  125. std::cout << "Loading VTK .vtu file " << fname << std::endl;
  126. vtkXMLUnstructuredGridReader* MeshReader = vtkXMLUnstructuredGridReader::New();
  127. MeshReader->SetFileName(fname.c_str());
  128. MeshReader->Update();
  129. //vtkGrid->DeepCopy( MeshReader->GetOutput() );
  130. //vtkGrid->ShallowCopy( MeshReader->GetOutput() );
  131. vtkGrid = MeshReader->GetOutput();
  132. MeshReader->Delete();
  133. return ;
  134. } // ----- end of method FEM4EllipticPDE::SetVTKGridFile -----
  135. //--------------------------------------------------------------------------------------
  136. // Class: FEM4EllipticPDE
  137. // Method: SetVTUGridFile
  138. //--------------------------------------------------------------------------------------
  139. void FEM4EllipticPDE::SetVTUGridFile ( char* fname ) {
  140. std::cout << "Loading VTK .vtu file " << fname << std::endl;
  141. vtkXMLUnstructuredGridReader* MeshReader = vtkXMLUnstructuredGridReader::New();
  142. MeshReader->SetFileName(fname);
  143. MeshReader->Update();
  144. //vtkGrid->DeepCopy( MeshReader->GetOutput() );
  145. //vtkGrid->ShallowCopy( MeshReader->GetOutput() );
  146. vtkGrid = MeshReader->GetOutput();
  147. MeshReader->Delete();
  148. return ;
  149. } // ----- end of method FEM4EllipticPDE::SetVTKGridFile -----
  150. vtkSmartPointer<vtkIdList> FEM4EllipticPDE::GetConnectedPoints(const int& id0) {
  151. vtkSmartPointer<vtkIdList> pointIds = vtkSmartPointer<vtkIdList>::New();
  152. vtkSmartPointer<vtkIdList> cellList = vtkSmartPointer<vtkIdList>::New();
  153. vtkGrid->GetPointCells(id0, cellList);
  154. for(int i=0;i<cellList->GetNumberOfIds(); ++i){
  155. vtkCell* cell = vtkGrid->GetCell(cellList->GetId(i));
  156. if(cell->GetNumberOfEdges() > 0){
  157. for(int j=0; j<cell->GetNumberOfEdges(); ++j){
  158. vtkCell* edge = cell->GetEdge(j);
  159. vtkIdList* edgePoints=edge->GetPointIds();
  160. if(edgePoints->GetId(0)==id0){
  161. pointIds->InsertUniqueId(edgePoints->GetId(1));
  162. } else if(edgePoints->GetId(1)==id0){
  163. pointIds->InsertUniqueId(edgePoints->GetId(0));
  164. }
  165. }
  166. }
  167. }
  168. return pointIds;
  169. }
  170. Real FEM4EllipticPDE::dist(Real r0[3], Real r1[3]) {
  171. Real rm0 = r1[0] - r0[0];
  172. Real rm1 = r1[1] - r0[1];
  173. Real rm2 = r1[2] - r0[2];
  174. return std::sqrt( rm0*rm0 + rm1*rm1 + rm2*rm2 );
  175. }
  176. Real FEM4EllipticPDE::dist(const Vector3r& r0, const Vector3r& r1) {
  177. Real rm0 = r1[0] - r0[0];
  178. Real rm1 = r1[1] - r0[1];
  179. Real rm2 = r1[2] - r0[2];
  180. return std::sqrt( rm0*rm0 + rm1*rm1 + rm2*rm2 );
  181. }
  182. //--------------------------------------------------------------------------------------
  183. // Class: FEM4EllipticPDE
  184. // Method: SetupDC
  185. //--------------------------------------------------------------------------------------
  186. void FEM4EllipticPDE::SetupDC ( DCSurvey* Survey, const int& ij ) {
  187. ////////////////////////////////////////////////////////////
  188. // Load vector g, solution vector u
  189. std::cout << "\nBuilding load vector (g)" << std::endl;
  190. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  191. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  192. int iia(0);
  193. Real jja(0);
  194. Survey->GetA( ij, iia, jja );
  195. //g(ii) = jj;
  196. int iib(0);
  197. Real jjb(0);
  198. Survey->GetB( ij, iib, jjb );
  199. //g(ii) = jj;
  200. /* 3D Phi */
  201. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  202. // Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  203. // for (int ii=0; ii<4; ++ii) {
  204. // double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  205. // C(ii, 0) = 1;
  206. // C(ii, 1) = pts[0];
  207. // C(ii, 2) = pts[1];
  208. // C(ii, 3) = pts[2];
  209. // }
  210. // Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  211. //
  212. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  213. int ID[4];
  214. ID[0] = Ids->GetId(0);
  215. ID[1] = Ids->GetId(1);
  216. ID[2] = Ids->GetId(2);
  217. ID[3] = Ids->GetId(3);
  218. //Real V = C.determinant(); // volume of tetrahedra
  219. Real sum(0);
  220. if (ID[0] == iia || ID[1] == iia || ID[2] == iia || ID[3] == iia ) {
  221. std::cout << "Caught A electrode, injecting " << iia << std::endl;
  222. //sum = 10;
  223. //g(ID[iia]) += jja;
  224. g(iia) += jja;
  225. }
  226. if (ID[0] == iib || ID[1] == iib || ID[2] == iib || ID[3] == iib) {
  227. //sum = -10;
  228. std::cout << "Caught B electrode, injecting " << iib << std::endl;
  229. //g(ID[iib]) += jjb;
  230. g(iib) += jjb;
  231. }
  232. //g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  233. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  234. }
  235. return ;
  236. } // ----- end of method FEM4EllipticPDE::SetupDC -----
  237. void FEM4EllipticPDE::Solve( const std::string& resfile ) {
  238. ConstructAMatrix();
  239. //SetupLineSourcePotential();
  240. if (vtkGrid->GetPointData()->GetScalars("G") ) {
  241. SetupSurfaceSourcePotential();
  242. } else if (vtkGrid->GetCellData()->GetScalars("G") ) {
  243. SetupPotential();
  244. }
  245. //ConstructLoadVector();
  246. std::cout << "\nSolving" << std::endl;
  247. std::cout << std::setw(5) << " " << std::setw(14) << "rows" << std::setw(14) << "cols" << std::endl;
  248. std::cout << std::setw(5) << " " << std::setw(14) << "--------" << std::setw(14) << "--------" << std::endl;
  249. std::cout << std::setw(5) << "A:" << std::setw(14) << A.rows() << std::setw(14) << A.cols() << std::endl;
  250. std::cout << std::setw(5) << "g:" << std::setw(14) << g.rows() << std::setw(14) << g.cols() << std::endl;
  251. ////////////////////////////////////////////////////////////
  252. // Solving:
  253. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  254. //VectorXr u = chol.solve(g);
  255. //#define LUSOLVE
  256. #ifdef LUSOLVE
  257. Eigen::SparseLU<Eigen::SparseMatrix<Real, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > solver;
  258. std::cout << "LU analyze pattern" << std::endl;
  259. solver.analyzePattern(A);
  260. std::cout << "LU factorizing" << std::endl;
  261. solver.factorize(A);
  262. std::cout << "LU solving" << std::endl;
  263. solver.factorize(A);
  264. VectorXr u = solver.solve(g);
  265. #endif
  266. #define CGSOLVE
  267. #ifdef CGSOLVE
  268. // TODO try IncompleteLUT preconditioner
  269. Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  270. //cg.setMaxIterations(3000);
  271. //cg.setTolerance(1e-28);
  272. VectorXr u = cg.solve(g);
  273. std::cout << "#iterations: " << cg.iterations() << std::endl;
  274. std::cout << "estimated error: " << cg.error() << std::endl;
  275. #endif
  276. vtkDoubleArray *gArray = vtkDoubleArray::New();
  277. vtkDoubleArray *uArray = vtkDoubleArray::New();
  278. uArray->SetNumberOfComponents(1);
  279. gArray->SetNumberOfComponents(1);
  280. for (int iu = 0; iu<u.size(); ++iu) {
  281. uArray->InsertTuple1(iu, u[iu]);
  282. gArray->InsertTuple1(iu, g[iu]);
  283. }
  284. uArray->SetName("u");
  285. gArray->SetName("g");
  286. vtkGrid->GetPointData()->AddArray(uArray);
  287. vtkGrid->GetPointData()->AddArray(gArray);
  288. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  289. Writer->SetInputData(vtkGrid);
  290. Writer->SetFileName(resfile.c_str());
  291. Writer->Write();
  292. Writer->Delete();
  293. gArray->Delete();
  294. uArray->Delete();
  295. }
  296. //--------------------------------------------------------------------------------------
  297. // Class: FEM4EllipticPDE
  298. // Method: ConstructAMatrix
  299. //--------------------------------------------------------------------------------------
  300. void FEM4EllipticPDE::ConstructAMatrix ( ) {
  301. /////////////////////////////////////////////////////////////////////////
  302. // Build stiffness matrix (A)
  303. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  304. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  305. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  306. //Eigen::SparseMatrix<Real>
  307. A.resize(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  308. std::vector< Eigen::Triplet<Real> > coeffs;
  309. if ( !vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet") ) {
  310. throw std::runtime_error("No HomogeneousDirichlet boundary conditions in input file.");
  311. }
  312. if ( !vtkGrid->GetCellData()->GetScalars("G") && !vtkGrid->GetPointData()->GetScalars("G") ) {
  313. throw std::runtime_error("No Cell or Point Data G");
  314. }
  315. bool GCell = false;
  316. if ( vtkGrid->GetCellData()->GetScalars("G") ) {
  317. GCell = true;
  318. }
  319. // Here we iterate over all of the cells
  320. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  321. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  322. // TODO, in production code we might not want to do this check here
  323. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  324. throw std::runtime_error("Non-tetrahedral mesh encountered!");
  325. }
  326. // construct coordinate matrix C
  327. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  328. for (int ii=0; ii<4; ++ii) {
  329. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  330. C(ii, 0) = 1;
  331. C(ii, 1) = pts[0] ;
  332. C(ii, 2) = pts[1] ;
  333. C(ii, 3) = pts[2] ;
  334. }
  335. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  336. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  337. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  338. int ID[4];
  339. ID[0] = Ids->GetId(0);
  340. ID[1] = Ids->GetId(1);
  341. ID[2] = Ids->GetId(2);
  342. ID[3] = Ids->GetId(3);
  343. Real sigma_bar(0);
  344. // TEST VOID IN SIGMA!! TODO DON"T KEEP THIS
  345. /*
  346. Real xc = C.col(1).array().mean();
  347. Real yc = C.col(2).array().mean();
  348. Real zc = C.col(3).array().mean();
  349. if ( xc >= 2.5 && xc <= 5. && yc>=2.5 && yc <= 5.) {
  350. sigma_bar = 0.;
  351. } else {
  352. sigma_bar = 1.;
  353. }
  354. */
  355. sigma_bar = 1.;
  356. for (int ii=0; ii<4; ++ii) {
  357. int bbi = vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0];
  358. if (bbi) {
  359. /* Dirichlet boundary */
  360. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[ii], 1));
  361. } else {
  362. for (int jj=0; jj<4; ++jj) {
  363. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  364. }
  365. }
  366. }
  367. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  368. }
  369. A.setFromTriplets(coeffs.begin(), coeffs.end());
  370. A.finalize();
  371. A.makeCompressed();
  372. }
  373. void FEM4EllipticPDE::SetupPotential() {
  374. ////////////////////////////////////////////////////////////
  375. // Load vector g
  376. std::cout << "\nBuilding load vector (g)" << std::endl;
  377. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  378. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " points" << std::endl;
  379. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  380. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  381. for (int ii=0; ii<4; ++ii) {
  382. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  383. C(ii, 0) = 1;
  384. C(ii, 1) = pts[0];
  385. C(ii, 2) = pts[1];
  386. C(ii, 3) = pts[2];
  387. }
  388. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  389. //Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  390. /* The indices */
  391. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  392. int ID[4];
  393. ID[0] = Ids->GetId(0);
  394. ID[1] = Ids->GetId(1);
  395. ID[2] = Ids->GetId(2);
  396. ID[3] = Ids->GetId(3);
  397. /* Fill the RHS vector with Dirichlet conditions or fuction value */
  398. for (int ii=0; ii<4; ++ii) {
  399. if (vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0]) {
  400. g(ID[ii]) += vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  401. } else {
  402. g(ID[ii]) += (V/4.)*(vtkGrid->GetCellData()->GetScalars("G")->GetTuple(ic)[0]); // Cell data
  403. /* for point data, determine if it is a point or line source */
  404. //g(ID[ii]) = (vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0]); // Point source
  405. //g(ID[ii]) = (vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0]); // Point source
  406. }
  407. }
  408. }
  409. }
  410. void FEM4EllipticPDE::SetupLineSourcePotential() {
  411. std::cerr << " FEM4EllipticPDE::SetupLineSourcePotential is not completed\n";
  412. exit(1);
  413. ////////////////////////////////////////////////////////////
  414. // Load vector g
  415. std::cout << "\nBuilding load vector (g)" << std::endl;
  416. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  417. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " points" << std::endl;
  418. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  419. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  420. for (int ii=0; ii<4; ++ii) {
  421. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  422. C(ii, 0) = 1;
  423. C(ii, 1) = pts[0];
  424. C(ii, 2) = pts[1];
  425. C(ii, 3) = pts[2];
  426. }
  427. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  428. //Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  429. /* The indices */
  430. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  431. int ID[4];
  432. ID[0] = Ids->GetId(0);
  433. ID[1] = Ids->GetId(1);
  434. ID[2] = Ids->GetId(2);
  435. ID[3] = Ids->GetId(3);
  436. /* Line source between nodes */
  437. for (int ii=0; ii<4; ++ii) {
  438. if (vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0]) {
  439. g(ID[ii]) += vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  440. } else {
  441. Real nv1 = vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0];
  442. for (int jj=ii+1; jj<4; ++jj) {
  443. Real nv2 = vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[jj])[0];
  444. if ( std::abs(nv1) > 1e-12 && std::abs(nv2) > 1e-12) {
  445. Real length = ( C.row(ii).tail<3>()-C.row(jj).tail<3>() ).norm();
  446. g(ID[ii]) += ((nv1+nv2)/2.)/(length);//*(V/4.);// * (vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0]); // Point source
  447. }
  448. }
  449. }
  450. }
  451. }
  452. }
  453. void FEM4EllipticPDE::SetupSurfaceSourcePotential() {
  454. ////////////////////////////////////////////////////////////
  455. // Load vector g
  456. std::cout << "\nBuilding load vector (g)" << std::endl;
  457. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  458. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  459. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  460. for (int ii=0; ii<4; ++ii) {
  461. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  462. C(ii, 0) = 1;
  463. C(ii, 1) = pts[0];
  464. C(ii, 2) = pts[1];
  465. C(ii, 3) = pts[2];
  466. }
  467. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  468. //Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  469. /* The indices */
  470. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  471. int ID[4];
  472. ID[0] = Ids->GetId(0);
  473. ID[1] = Ids->GetId(1);
  474. ID[2] = Ids->GetId(2);
  475. ID[3] = Ids->GetId(3);
  476. // Apply Dirichlet conditions
  477. for (int ii=0; ii<4; ++ii) {
  478. //if (vtkGrid->GetPointData()->GetScalars("InHomogeneousDirichlet")->GetTuple(ID[ii])[0]) {
  479. // g(ID[ii]) += vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  480. //}
  481. }
  482. /* the 4 faces of the tetrahedra
  483. ID[0] ID[1] ID[2]
  484. ID[0] ID[1] ID[3]
  485. ID[0] ID[2] ID[3]
  486. ID[1] ID[2] ID[3]
  487. */
  488. Real i4pi = .5;
  489. /* surfaces of tetrahedra */
  490. // Face 0, ID 0,1,2
  491. Eigen::Matrix<Real, 3, 3> CC = Eigen::Matrix<Real, 3, 3>::Ones() ;
  492. if ( std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])) > 1e-12 &&
  493. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])) > 1e-12 &&
  494. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])) > 1e-12 ) {
  495. CC.col(1) = C.row(0).tail<3>() - C.row(1).tail<3>();
  496. CC.col(2) = C.row(0).tail<3>() - C.row(2).tail<3>();
  497. Real TA = .5*CC.col(1).cross(CC.col(2)).norm();
  498. g(ID[0]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])*TA/3. * i4pi ;
  499. g(ID[1]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])*TA/3. * i4pi ;
  500. g(ID[2]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])*TA/3. * i4pi ;
  501. }
  502. // Face 1, ID 0,1,3
  503. if ( std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])) > 1e-12 &&
  504. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])) > 1e-12 &&
  505. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])) > 1e-12 ) {
  506. CC.col(1) = C.row(0).tail<3>() - C.row(1).tail<3>();
  507. CC.col(2) = C.row(0).tail<3>() - C.row(3).tail<3>();
  508. Real TA = .5*CC.col(1).cross(CC.col(2)).norm();
  509. g(ID[0]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])*TA/3. * i4pi ;
  510. g(ID[1]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])*TA/3. * i4pi ;
  511. g(ID[3]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])*TA/3. * i4pi ;
  512. }
  513. // Face 2, ID 0,2,3
  514. if ( std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])) > 1e-12 &&
  515. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])) > 1e-12 &&
  516. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])) > 1e-12 ) {
  517. CC.col(1) = C.row(0).tail<3>() - C.row(2).tail<3>();
  518. CC.col(2) = C.row(0).tail<3>() - C.row(3).tail<3>();
  519. Real TA = .5*CC.col(1).cross(CC.col(2)).norm();
  520. g(ID[0]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0])*TA/3. * i4pi ;
  521. g(ID[2]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])*TA/3. * i4pi ;
  522. g(ID[3]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])*TA/3. * i4pi ;
  523. }
  524. // Face 3, ID 1,2,3
  525. if ( std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])) > 1e-12 &&
  526. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])) > 1e-12 &&
  527. std::abs(vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])) > 1e-12 ) {
  528. CC.col(1) = C.row(1).tail<3>() - C.row(2).tail<3>();
  529. CC.col(2) = C.row(1).tail<3>() - C.row(3).tail<3>();
  530. Real TA = .5*CC.col(1).cross(CC.col(2)).norm();
  531. g(ID[1]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1])*TA/3. * i4pi ;
  532. g(ID[2]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2])*TA/3. * i4pi ;
  533. g(ID[3]) += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3])*TA/3. * i4pi ;
  534. }
  535. }
  536. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " points" << std::endl;
  537. }
  538. #if 0
  539. void FEM4EllipticPDE::SolveOLD2(const std::string& fname) {
  540. Real r0[3];
  541. Real r1[3];
  542. /////////////////////////////////////////////////////////////////////////
  543. // Surface filter, to determine if points are on boundary, and need
  544. // boundary conditions applied
  545. vtkDataSetSurfaceFilter* Surface = vtkDataSetSurfaceFilter::New();
  546. Surface->SetInputData(vtkGrid);
  547. Surface->PassThroughPointIdsOn( );
  548. Surface->Update();
  549. vtkIdTypeArray* BdryIds = static_cast<vtkIdTypeArray*>
  550. (Surface->GetOutput()->GetPointData()->GetScalars("vtkOriginalPointIds"));
  551. // Expensive search for whether or not point is on boundary. O(n) cost.
  552. VectorXi bndryCnt = VectorXi::Zero(vtkGrid->GetNumberOfPoints());
  553. for (int isp=0; isp < Surface->GetOutput()->GetNumberOfPoints(); ++isp) {
  554. //double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  555. // x \in -14.5 to 14.5
  556. // y \in 0 to 30
  557. bndryCnt(BdryIds->GetTuple1(isp)) += 1;
  558. }
  559. /////////////////////////////////////////////////////////////////////////
  560. // Build stiffness matrix (A)
  561. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  562. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  563. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  564. Eigen::SparseMatrix<Real> A(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  565. std::vector< Eigen::Triplet<Real> > coeffs;
  566. // Here we iterate over all of the cells
  567. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  568. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  569. // TODO, in production code we might not want to do this check here
  570. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  571. std::cout << "DOOM FEM4EllipticPDE encountered non-tetrahedral mesh\n";
  572. std::cout << "Number of points in cell " << vtkGrid->GetCell(ic)->GetNumberOfPoints() << std::endl ;
  573. exit(1);
  574. }
  575. // construct coordinate matrix C
  576. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  577. for (int ii=0; ii<4; ++ii) {
  578. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  579. C(ii, 0) = 1;
  580. C(ii, 1) = pts[0] ;
  581. C(ii, 2) = pts[1] ;
  582. C(ii, 3) = pts[2] ;
  583. }
  584. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  585. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  586. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  587. int ID[4];
  588. ID[0] = Ids->GetId(0);
  589. ID[1] = Ids->GetId(1);
  590. ID[2] = Ids->GetId(2);
  591. ID[3] = Ids->GetId(3);
  592. Real sum(0);
  593. Real sigma_bar = vtkGrid->GetCellData()->GetScalars()->GetTuple1(ic);
  594. for (int ii=0; ii<4; ++ii) {
  595. for (int jj=0; jj<4; ++jj) {
  596. if (jj == ii) {
  597. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  598. // solve for the boundaries? Is one better? This seems to work, which is nice.
  599. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ii] ); // + sum;
  600. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ii])[0];
  601. //std::cout << "bb " << bb << std::endl;
  602. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  603. 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 ) );
  604. } else {
  605. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  606. }
  607. // Stiffness matrix no longer contains boundary conditions...
  608. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  609. }
  610. }
  611. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  612. }
  613. A.setFromTriplets(coeffs.begin(), coeffs.end());
  614. //A.makeCompressed();
  615. ////////////////////////////////////////////////////////////
  616. // Load vector g, solution vector u
  617. std::cout << "\nBuilding load vector (g)" << std::endl;
  618. VectorXr g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  619. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  620. // If the G function has been evaluated at each *node*
  621. // --> but still needs to be integrated at the surfaces
  622. // Aha, requires that there is in fact a pointdata memeber // BUG TODO BUG!!!
  623. std::cout << "Point Data ptr " << vtkGrid->GetPointData() << std::endl;
  624. //if ( vtkGrid->GetPointData() != NULL && std::string( vtkGrid->GetPointData()->GetScalars()->GetName() ).compare( std::string("G") ) == 0 ) {
  625. bool pe(false);
  626. bool ne(false);
  627. if ( true ) {
  628. std::cout << "\nUsing G from file" << std::endl;
  629. /* 3D Phi */
  630. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  631. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  632. for (int ii=0; ii<4; ++ii) {
  633. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  634. C(ii, 0) = 1;
  635. C(ii, 1) = pts[0];
  636. C(ii, 2) = pts[1];
  637. C(ii, 3) = pts[2];
  638. }
  639. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  640. //Real V = C.determinant(); // volume of tetrahedra
  641. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  642. int ID[4];
  643. ID[0] = Ids->GetId(0);
  644. ID[1] = Ids->GetId(1);
  645. ID[2] = Ids->GetId(2);
  646. ID[3] = Ids->GetId(3);
  647. /* bad news bears for magnet */
  648. double* pt = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(0);
  649. Real sum(0);
  650. /*
  651. if (!pe) {
  652. if (std::abs(pt[0]) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  653. sum = 1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  654. pe = true;
  655. }
  656. }*/
  657. if (ID[0] == 26) {
  658. sum = 10;
  659. }
  660. if (ID[0] == 30) {
  661. sum = -10;
  662. }
  663. /*
  664. if (!ne) {
  665. if (std::abs(pt[0]+1.) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  666. sum = -1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  667. std::cout << "Negative Electroce\n";
  668. ne = true;
  669. }
  670. }
  671. */
  672. //for (int ii=0; ii<4; ++ii) {
  673. //g(ID[ii]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0] ) ;
  674. //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0]) > 1e-3 )
  675. //}
  676. // TODO check Load Vector...
  677. g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  678. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  679. }
  680. /*
  681. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  682. vtkGrid->GetPoint(ic, r0);
  683. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  684. double g0 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ;
  685. //std::cout << "num conn " << connectedVertices->GetNumberOfIds() << std::endl;
  686. if ( std::abs(g0) > 1e-3 ) {
  687. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  688. int ii = connectedVertices->GetId(i);
  689. vtkGrid->GetPoint(ii, r1);
  690. double g1 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ii)[0] ;
  691. //g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  692. if ( std::abs(g1) > 1e-3 ) {
  693. g(ic) += CompositeSimpsons2(g1, g0, r1, r0, 1000);
  694. }
  695. //g(ic) += CompositeSimpsons2(g0, r1, r0, 8);
  696. //if ( std::abs(g1) > 1e-3 ) {
  697. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 8);
  698. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 100); // / (2*dist(r0,r1)) ;
  699. // g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  700. //g(ic) += CompositeSimpsons2(g0, r0, r1, 8);
  701. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  702. //} //else {
  703. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  704. //}
  705. }
  706. }
  707. //g(ic) = 2.* vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ; // Why 2?
  708. //std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  709. }
  710. */
  711. } else if (vtkG) { // VTK implicit function, proceed with care
  712. std::cout << "\nUsing implicit file from file" << std::endl;
  713. // OpenMP right here
  714. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  715. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  716. //vtkGrid->GetPoint(ic, r0);
  717. //g(ic) += vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  718. // std::cout << vtkG->FunctionValue(r0[0], r0[1], r0[2]) << std::endl;
  719. //g(ic) += vtkGrid->GetPointData()->GetScalars()->GetTuple1(ic);// FunctionValue(r0[0], r0[1], r0[2]) ;
  720. /*
  721. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  722. int ii = connectedVertices->GetId(i);
  723. vtkGrid->GetPoint(ii, r1);
  724. g(ic) += CompositeSimpsons2(vtkG, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  725. }
  726. */
  727. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  728. }
  729. } else if (gFcn3) {
  730. std::cout << "\nUsing gFcn3" << std::endl;
  731. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  732. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  733. vtkGrid->GetPoint(ic, r0);
  734. // TODO, test OMP sum reduction here. Is vtkGrid->GetPoint thread safe?
  735. //Real sum(0.);
  736. //#ifdef LEMMAUSEOMP
  737. //#pragma omp parallel for reduction(+:sum)
  738. //#endif
  739. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  740. int ii = connectedVertices->GetId(i);
  741. vtkGrid->GetPoint(ii, r1);
  742. g(ic) += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  743. //sum += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  744. }
  745. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  746. //g(ic) = sum;
  747. }
  748. } else {
  749. std::cout << "No source specified\n";
  750. exit(EXIT_FAILURE);
  751. }
  752. // std::cout << g << std::endl;
  753. //g(85) = 1;
  754. std::cout << "\nSolving" << std::endl;
  755. ////////////////////////////////////////////////////////////
  756. // Solving:
  757. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  758. //VectorXr u = chol.solve(g);
  759. //Eigen::SparseLU<Eigen::SparseMatrix<Real, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > solver;
  760. //solver.analyzePattern(A);
  761. //solver.factorize(A);
  762. //VectorXr u = solver.solve(g);
  763. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  764. Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  765. //Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  766. cg.setMaxIterations(3000);
  767. //cg.compute(A);
  768. //std::cout << "Computed " << std::endl;
  769. VectorXr u = cg.solve(g);
  770. std::cout << "#iterations: " << cg.iterations() << std::endl;
  771. std::cout << "estimated error: " << cg.error() << std::endl;
  772. vtkDoubleArray *gArray = vtkDoubleArray::New();
  773. vtkDoubleArray *uArray = vtkDoubleArray::New();
  774. uArray->SetNumberOfComponents(1);
  775. gArray->SetNumberOfComponents(1);
  776. for (int iu = 0; iu<u.size(); ++iu) {
  777. uArray->InsertTuple1(iu, u[iu]);
  778. gArray->InsertTuple1(iu, g[iu]);
  779. }
  780. uArray->SetName("u");
  781. gArray->SetName("g");
  782. vtkGrid->GetPointData()->AddArray(uArray);
  783. vtkGrid->GetPointData()->AddArray(gArray);
  784. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  785. Writer->SetInputData(vtkGrid);
  786. Writer->SetFileName(fname.c_str());
  787. Writer->Write();
  788. Writer->Delete();
  789. Surface->Delete();
  790. gArray->Delete();
  791. uArray->Delete();
  792. }
  793. #endif
  794. // Uses simpon's rule to perform a definite integral of a
  795. // function (passed as a pointer). The integration occurs from
  796. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  797. Real FEM4EllipticPDE::CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int n) {
  798. Vector3r R0(r0[0], r0[1], r0[2]);
  799. Vector3r R1(r1[0], r1[1], r1[2]);
  800. // force n to be even
  801. assert(n > 0);
  802. //n += (n % 2);
  803. Real h = dist(r0, r1) / (Real)(n) ;
  804. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  805. Vector3r dr = (R1 - R0).array() / Real(n);
  806. Vector3r rx;
  807. rx.array() = R0.array() + dr.array();
  808. for (int i=1; i<n; i+=2) {
  809. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]);
  810. rx += 2.*dr;
  811. }
  812. rx.array() = R0.array() + 2*dr.array();
  813. for (int i=2; i<n-1; i+=2) {
  814. S += 2.*f->FunctionValue(rx[0], rx[1], rx[2]);
  815. rx += 2.*dr;
  816. }
  817. return h * S / 3.;
  818. }
  819. // Uses simpon's rule to perform a definite integral of a
  820. // function (passed as a pointer). The integration occurs from
  821. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  822. // This is just here as a convenience
  823. Real FEM4EllipticPDE::CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int n) {
  824. return dist(r0,r1)*f;
  825. /*
  826. Vector3r R0(r0[0], r0[1], r0[2]);
  827. Vector3r R1(r1[0], r1[1], r1[2]);
  828. // force n to be even
  829. assert(n > 0);
  830. //n += (n % 2);
  831. Real h = dist(r0, r1) / (Real)(n) ;
  832. Real S = f + f;
  833. Vector3r dr = (R1 - R0).array() / Real(n);
  834. //Vector3r rx;
  835. //rx.array() = R0.array() + dr.array();
  836. for (int i=1; i<n; i+=2) {
  837. S += 4. * f;
  838. //rx += 2.*dr;
  839. }
  840. //rx.array() = R0.array() + 2*dr.array();
  841. for (int i=2; i<n-1; i+=2) {
  842. S += 2. * f;
  843. //rx += 2.*dr;
  844. }
  845. return h * S / 3.;
  846. */
  847. }
  848. /*
  849. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  850. * test function owned by the FEM implimentaion.
  851. */
  852. Real FEM4EllipticPDE::CompositeSimpsons2(vtkImplicitFunction* f,
  853. Real r0[3], Real r1[3], int n) {
  854. Vector3r R0(r0[0], r0[1], r0[2]);
  855. Vector3r R1(r1[0], r1[1], r1[2]);
  856. // force n to be even
  857. assert(n > 0);
  858. //n += (n % 2);
  859. Real h = dist(r0, r1) / (Real)(n) ;
  860. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  861. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  862. //Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  863. Vector3r dr = (R1 - R0).array() / Real(n);
  864. Vector3r rx;
  865. rx.array() = R0.array() + dr.array();
  866. for (int i=1; i<n; i+=2) {
  867. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  868. rx += 2.*dr;
  869. }
  870. rx.array() = R0.array() + 2*dr.array();
  871. for (int i=2; i<n-1; i+=2) {
  872. S += 2. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  873. rx += 2.*dr;
  874. }
  875. return h * S / 3.;
  876. }
  877. /*
  878. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  879. * test function owned by the FEM implimentaion.
  880. */
  881. Real FEM4EllipticPDE::CompositeSimpsons2( Real (*f)(const Real&, const Real&, const Real&),
  882. Real r0[3], Real r1[3], int n) {
  883. Vector3r R0(r0[0], r0[1], r0[2]);
  884. Vector3r R1(r1[0], r1[1], r1[2]);
  885. // force n to be even
  886. assert(n > 0);
  887. //n += (n % 2);
  888. Real h = dist(r0, r1) / (Real)(n) ;
  889. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  890. //Real S = f(r0[0], r0[1], r0[2])*Hat(R0, R0, R1) + f(r1[0], r1[1], r1[2])*Hat(R1, R0, R1);
  891. Real S = f(r0[0], r0[1], r0[2]) + f(r1[0], r1[1], r1[2]);
  892. Vector3r dr = (R1 - R0).array() / Real(n);
  893. Vector3r rx;
  894. rx.array() = R0.array() + dr.array();
  895. for (int i=1; i<n; i+=2) {
  896. S += 4. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  897. rx += 2.*dr;
  898. }
  899. rx.array() = R0.array() + 2*dr.array();
  900. for (int i=2; i<n-1; i+=2) {
  901. S += 2. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  902. rx += 2.*dr;
  903. }
  904. return h * S / 3.;
  905. }
  906. /*
  907. * Performs numerical integration of two functions, one is constant valued f, the other is the FEM
  908. * test function owned by the FEM implimentaion.
  909. */
  910. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int n) {
  911. Vector3r R0(r0[0], r0[1], r0[2]);
  912. Vector3r R1(r1[0], r1[1], r1[2]);
  913. // force n to be even
  914. assert(n > 0);
  915. //n += (n % 2);
  916. Real h = dist(r0, r1) / (Real)(n) ;
  917. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  918. Real S = 2*f; //*Hat(R0, R0, R1) + f*Hat(R1, R0, R1);
  919. Vector3r dr = (R1 - R0).array() / Real(n);
  920. Vector3r rx;
  921. rx.array() = R0.array() + dr.array();
  922. for (int i=1; i<n; i+=2) {
  923. S += 4. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  924. rx += 2.*dr;
  925. }
  926. rx.array() = R0.array() + 2*dr.array();
  927. for (int i=2; i<n-1; i+=2) {
  928. S += 2. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  929. rx += 2.*dr;
  930. }
  931. return h * S / 3.;
  932. }
  933. /*
  934. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  935. * test function owned by the FEM implimentaion.
  936. */
  937. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  938. Vector3r R0(r0[0], r0[1], r0[2]);
  939. Vector3r R1(r1[0], r1[1], r1[2]);
  940. // force n to be even
  941. assert(n > 0);
  942. //n += (n % 2);
  943. Real h = dist(r0, r1) / (Real)(n) ;
  944. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  945. // NO! We are looking at 1/2 hat often!!! So S = f1!
  946. Real S = f1; //f0*Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  947. Vector3r dr = (R1 - R0).array() / Real(n);
  948. // linear interpolate function
  949. //Vector3r rx;
  950. //rx.array() = R0.array() + dr.array();
  951. for (int i=1; i<n; i+=2) {
  952. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  953. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0) ;
  954. }
  955. //rx.array() = R0.array() + 2*dr.array();
  956. for (int i=2; i<n-1; i+=2) {
  957. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  958. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0);
  959. }
  960. return h * S / 3.;
  961. }
  962. /*
  963. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  964. * test function owned by the FEM implimentaion.
  965. */
  966. Real FEM4EllipticPDE::CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  967. Vector3r R0(r0[0], r0[1], r0[2]);
  968. Vector3r R1(r1[0], r1[1], r1[2]);
  969. // force n to be even
  970. assert(n > 0);
  971. //n += (n % 2);
  972. Real h = dist(r0, r1) / (Real)(n) ;
  973. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  974. // NO! We are looking at 1/2 hat often!!! So S = f1!
  975. Real S = f0+f1; //Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  976. Vector3r dr = (R1 - R0).array() / Real(n);
  977. // linear interpolate function
  978. //Vector3r rx;
  979. //rx.array() = R0.array() + dr.array();
  980. for (int i=1; i<n; i+=2) {
  981. double fx = 1;// f0 + (f1 - f0) * ((i*h)/(h*n));
  982. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1) * Hat(R1.array() + i*dr.array(), R1, R0) ;
  983. }
  984. //rx.array() = R0.array() + 2*dr.array();
  985. for (int i=2; i<n-1; i+=2) {
  986. double fx = 1; // f0 + (f1 - f0) * ((i*h)/(h*n));
  987. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1)* Hat(R1.array() + i*dr.array(), R1, R0);
  988. }
  989. return h * S / 3.;
  990. }
  991. //--------------------------------------------------------------------------------------
  992. // Class: FEM4EllipticPDE
  993. // Method: Hat
  994. //--------------------------------------------------------------------------------------
  995. Real FEM4EllipticPDE::Hat ( const Vector3r& r, const Vector3r& r0, const Vector3r& r1 ) {
  996. //return (r-r0).norm() / (r1-r0).norm() ;
  997. return dist(r, r0) / dist(r1, r0) ;
  998. } // ----- end of method FEM4EllipticPDE::Hat -----
  999. } // ----- end of Lemma name -----