Galerkin FEM for elliptic PDEs
Você não pode selecionar mais de 25 tópicos Os tópicos devem começar com uma letra ou um número, podem incluir traços ('-') e podem ter até 35 caracteres.

FEM4EllipticPDE.cpp 38KB

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