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

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