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

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