Galerkin FEM for elliptic PDEs
Nelze vybrat více než 25 témat Téma musí začínat písmenem nebo číslem, může obsahovat pomlčky („-“) a může být dlouhé až 35 znaků.

FEM4EllipticPDE.cpp 41KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014
  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. SetupPotential();
  167. //ConstructLoadVector();
  168. std::cout << "\nSolving" << std::endl;
  169. std::cout << " rows\tcols\n";
  170. std::cout << "A: " << A.rows() << "\t" << A.cols() << std::endl;
  171. std::cout << "g: " << g.rows() << "\t" << g.cols() << std::endl;
  172. ////////////////////////////////////////////////////////////
  173. // Solving:
  174. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  175. //VectorXr u = chol.solve(g);
  176. //Eigen::SparseLU<Eigen::SparseMatrix<Real, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > solver;
  177. //solver.analyzePattern(A);
  178. //solver.factorize(A);
  179. //VectorXr u = solver.solve(g);
  180. //#ifdef CGSOLVE
  181. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  182. //Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  183. Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  184. //cg.setMaxIterations(3000);
  185. //cg.setTolerance(1e-28);
  186. VectorXr u = cg.solve(g);
  187. std::cout << "#iterations: " << cg.iterations() << std::endl;
  188. std::cout << "estimated error: " << cg.error() << std::endl;
  189. //#endif
  190. vtkDoubleArray *gArray = vtkDoubleArray::New();
  191. vtkDoubleArray *uArray = vtkDoubleArray::New();
  192. uArray->SetNumberOfComponents(1);
  193. gArray->SetNumberOfComponents(1);
  194. for (int iu = 0; iu<u.size(); ++iu) {
  195. uArray->InsertTuple1(iu, u[iu]);
  196. gArray->InsertTuple1(iu, g[iu]);
  197. }
  198. uArray->SetName("u");
  199. gArray->SetName("g");
  200. vtkGrid->GetPointData()->AddArray(uArray);
  201. vtkGrid->GetPointData()->AddArray(gArray);
  202. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  203. Writer->SetInputData(vtkGrid);
  204. Writer->SetFileName(resfile.c_str());
  205. Writer->Write();
  206. Writer->Delete();
  207. gArray->Delete();
  208. uArray->Delete();
  209. }
  210. //--------------------------------------------------------------------------------------
  211. // Class: FEM4EllipticPDE
  212. // Method: ConstructAMatrix
  213. //--------------------------------------------------------------------------------------
  214. void FEM4EllipticPDE::ConstructAMatrix ( ) {
  215. /////////////////////////////////////////////////////////////////////////
  216. // Build stiffness matrix (A)
  217. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  218. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  219. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  220. //Eigen::SparseMatrix<Real>
  221. A.resize(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  222. std::vector< Eigen::Triplet<Real> > coeffs;
  223. if ( !vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet") ) {
  224. throw std::runtime_error("No HomogeneousDirichlet boundary conditions in input file.");
  225. }
  226. if ( !vtkGrid->GetCellData()->GetScalars("G") && !vtkGrid->GetPointData()->GetScalars("G") ) {
  227. throw std::runtime_error("No Cell or Point Data G");
  228. }
  229. bool GCell = false;
  230. if ( vtkGrid->GetCellData()->GetScalars("G") ) {
  231. GCell = true;
  232. }
  233. // Here we iterate over all of the cells
  234. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  235. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  236. // TODO, in production code we might not want to do this check here
  237. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  238. throw std::runtime_error("Non-tetrahedral mesh encountered!");
  239. }
  240. // construct coordinate matrix C
  241. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  242. for (int ii=0; ii<4; ++ii) {
  243. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  244. C(ii, 0) = 1;
  245. C(ii, 1) = pts[0] ;
  246. C(ii, 2) = pts[1] ;
  247. C(ii, 3) = pts[2] ;
  248. }
  249. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  250. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  251. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  252. int ID[4];
  253. ID[0] = Ids->GetId(0);
  254. ID[1] = Ids->GetId(1);
  255. ID[2] = Ids->GetId(2);
  256. ID[3] = Ids->GetId(3);
  257. Real sigma_bar(0);
  258. /*
  259. if (GCell) {
  260. sigma_bar = vtkGrid->GetCellData()->GetScalars("G")->GetTuple1(ic);
  261. } else {
  262. sigma_bar = vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[0]);
  263. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[1]);
  264. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[2]);
  265. sigma_bar += vtkGrid->GetPointData()->GetScalars("G")->GetTuple1(ID[3]);
  266. sigma_bar /= 4.;
  267. }
  268. */
  269. // TEST VOID IN SIGMA!! TODO DON"T KEEP THIS
  270. Real xc = C.col(1).array().mean();
  271. Real yc = C.col(2).array().mean();
  272. Real zc = C.col(3).array().mean();
  273. if ( xc >= 2.5 && xc <= 5. && yc>=2.5 && yc <= 5.) {
  274. sigma_bar = 0.;
  275. } else {
  276. sigma_bar = 1.;
  277. }
  278. sigma_bar = 1.;
  279. for (int ii=0; ii<4; ++ii) {
  280. int bbi = vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0];
  281. if (bbi) {
  282. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[ii], 1));
  283. } else {
  284. for (int jj=0; jj<4; ++jj) {
  285. /* homogeneous Dirichlet boundary */
  286. //if (jj==ii && C(ii,3)==13.5) {
  287. //if (jj==ii && ((bb-1.)<1e-8) ) {
  288. // Apply Homogeneous Dirichlet Boundary conditions
  289. // Real bb = vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0];
  290. // Real bdry = (1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  291. //Real bdry = GradPhi.col(ii).tail<3>().dot(GradPhi.col(ii).tail<3>())*BndrySigma*bb; // + sum;
  292. // 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 ) );
  293. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], 1));
  294. //break;
  295. //} else {
  296. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  297. //}
  298. /* */
  299. /* Inhomogeneous Dirichlet */
  300. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>())*V*sigma_bar ));
  301. // Stiffness matrix no longer contains boundary conditions...
  302. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  303. }
  304. }
  305. }
  306. //std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  307. }
  308. A.setFromTriplets(coeffs.begin(), coeffs.end());
  309. //std::cout << "A\n" << A << std::endl;
  310. }
  311. void FEM4EllipticPDE::SetupPotential() {
  312. ////////////////////////////////////////////////////////////
  313. // Load vector g
  314. std::cout << "\nBuilding load vector (g)" << std::endl;
  315. g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  316. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " points" << std::endl;
  317. VectorXr DB = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  318. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  319. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  320. for (int ii=0; ii<4; ++ii) {
  321. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  322. C(ii, 0) = 1;
  323. C(ii, 1) = pts[0];
  324. C(ii, 2) = pts[1];
  325. C(ii, 3) = pts[2];
  326. }
  327. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  328. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  329. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  330. int ID[4];
  331. ID[0] = Ids->GetId(0);
  332. ID[1] = Ids->GetId(1);
  333. ID[2] = Ids->GetId(2);
  334. ID[3] = Ids->GetId(3);
  335. /*
  336. Real avg(0);
  337. Real GG[4];
  338. for (int ii=0; ii<4; ++ii) {
  339. GG[ii] = vtkGrid->GetPointData()->GetScalars("G")->GetTuple(ID[ii])[0];
  340. //avg /= 4.;
  341. }
  342. if ( std::abs( (GG[0]+GG[1]+GG[2]+GG[3])/4. - GG[0]) < 1e-5) {
  343. avg = GG[0];
  344. }
  345. */
  346. /*
  347. VectorXr W = VectorXr::Zero(4);
  348. for (int ii=0; ii<4; ++ii) {
  349. W[ii] = vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0] *
  350. vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  351. DB[ID[ii]] = vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0] *
  352. vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  353. }
  354. //auto G = GradPhi.block<3,4>(1,0).transpose()*GradPhi.block<3,4>(1,0)*W;
  355. VectorXr G = GradPhi.block<3,4>(1,0).transpose()*GradPhi.block<3,4>(1,0)*W;
  356. Real sigma_bar(1.);
  357. */
  358. for (int ii=0; ii<4; ++ii) {
  359. // avg += vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  360. // //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0]) > 1e-3 )
  361. //}
  362. //TODO this seems wrong!
  363. //avg /= 4.;
  364. // TODO test code remove
  365. if (vtkGrid->GetPointData()->GetScalars("HomogeneousDirichlet")->GetTuple(ID[ii])[0]) {
  366. g(ID[ii]) += vtkGrid->GetPointData()->GetScalars("analytic_phi")->GetTuple(ID[ii])[0];
  367. } else {
  368. g(ID[ii]) += (PI*V)*(vtkGrid->GetCellData()->GetScalars("G")->GetTuple(ic)[0]); // Why 3.0??
  369. }
  370. //if (std::abs(C(ii,3)-13.5) > 1e-5) {
  371. // g(ID[ii]) -= G[ii]*(V)*sigma_bar;
  372. //}
  373. //g(ID[ii]) += V/4*avg;
  374. //g(ID[ii]) += 6.67 *(V/4.) * avg;
  375. }
  376. //g(ID[0]) += (V/4.) * avg;
  377. }
  378. //g -= A*DB;
  379. }
  380. void FEM4EllipticPDE::SolveOLD(const std::string& fname) {
  381. Real r0[3];
  382. Real r1[3];
  383. /////////////////////////////////////////////////////////////////////////
  384. // Surface filter, to determine if points are on boundary, and need
  385. // boundary conditions applied
  386. vtkDataSetSurfaceFilter* Surface = vtkDataSetSurfaceFilter::New();
  387. Surface->SetInputData(vtkGrid);
  388. Surface->PassThroughPointIdsOn( );
  389. Surface->Update();
  390. vtkIdTypeArray* BdryIds = static_cast<vtkIdTypeArray*>
  391. (Surface->GetOutput()->GetPointData()->GetScalars("vtkOriginalPointIds"));
  392. // Expensive search for whether or not point is on boundary. O(n) cost.
  393. VectorXi bndryCnt = VectorXi::Zero(vtkGrid->GetNumberOfPoints());
  394. for (int isp=0; isp < Surface->GetOutput()->GetNumberOfPoints(); ++isp) {
  395. //double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  396. // x \in -14.5 to 14.5
  397. // y \in 0 to 30
  398. bndryCnt(BdryIds->GetTuple1(isp)) += 1;
  399. }
  400. /////////////////////////////////////////////////////////////////////////
  401. // Build stiffness matrix (A)
  402. std::cout << "Building Stiffness Matrix (A) " << std::endl;
  403. std::cout << "\tMesh has " << vtkGrid->GetNumberOfCells() << " cells " << std::endl;
  404. std::cout << "\tMesh has " << vtkGrid->GetNumberOfPoints() << " points " << std::endl;
  405. Eigen::SparseMatrix<Real> A(vtkGrid->GetNumberOfPoints(), vtkGrid->GetNumberOfPoints());
  406. std::vector< Eigen::Triplet<Real> > coeffs;
  407. // Here we iterate over all of the cells
  408. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  409. assert ( vtkGrid->GetCell(ic)->GetNumberOfPoints() == 4 );
  410. // TODO, in production code we might not want to do this check here
  411. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  412. std::cout << "DOOM FEM4EllipticPDE encountered non-tetrahedral mesh\n";
  413. std::cout << "Number of points in cell " << vtkGrid->GetCell(ic)->GetNumberOfPoints() << std::endl ;
  414. exit(1);
  415. }
  416. // construct coordinate matrix C
  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. Eigen::Matrix<Real, 4, 4> GradPhi = C.inverse(); // nabla \phi
  426. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  427. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  428. int ID[4];
  429. ID[0] = Ids->GetId(0);
  430. ID[1] = Ids->GetId(1);
  431. ID[2] = Ids->GetId(2);
  432. ID[3] = Ids->GetId(3);
  433. Real sum(0);
  434. Real sigma_bar = vtkGrid->GetCellData()->GetScalars()->GetTuple1(ic);
  435. for (int ii=0; ii<4; ++ii) {
  436. for (int jj=0; jj<4; ++jj) {
  437. if (jj == ii) {
  438. // I apply boundary to Stiffness matrix, it's common to take the other approach and apply to the load vector and then
  439. // solve for the boundaries? Is one better? This seems to work, which is nice.
  440. //Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bndryCnt( ID[ii] ); // + sum;
  441. Real bb = vtkGrid->GetPointData()->GetScalars("vtkValidPointMask")->GetTuple(ID[ii])[0];
  442. //std::cout << "bb " << bb << std::endl;
  443. Real bdry = V*(1./(BndryH*BndryH))*BndrySigma*bb; // + sum;
  444. 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 ) );
  445. } else {
  446. coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  447. }
  448. // Stiffness matrix no longer contains boundary conditions...
  449. //coeffs.push_back( Eigen::Triplet<Real> ( ID[ii], ID[jj], GradPhi.col(ii).tail<3>().dot(GradPhi.col(jj).tail<3>() ) * V * sigma_bar ) );
  450. }
  451. }
  452. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  453. }
  454. A.setFromTriplets(coeffs.begin(), coeffs.end());
  455. //A.makeCompressed();
  456. ////////////////////////////////////////////////////////////
  457. // Load vector g, solution vector u
  458. std::cout << "\nBuilding load vector (g)" << std::endl;
  459. VectorXr g = VectorXr::Zero(vtkGrid->GetNumberOfPoints());
  460. std::cout << "made g with " << vtkGrid->GetNumberOfPoints() << " pnts" << std::endl;
  461. // If the G function has been evaluated at each *node*
  462. // --> but still needs to be integrated at the surfaces
  463. // Aha, requires that there is in fact a pointdata memeber // BUG TODO BUG!!!
  464. std::cout << "Point Data ptr " << vtkGrid->GetPointData() << std::endl;
  465. //if ( vtkGrid->GetPointData() != NULL && std::string( vtkGrid->GetPointData()->GetScalars()->GetName() ).compare( std::string("G") ) == 0 ) {
  466. bool pe(false);
  467. bool ne(false);
  468. if ( true ) {
  469. std::cout << "\nUsing G from file" << std::endl;
  470. /* 3D Phi */
  471. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  472. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  473. for (int ii=0; ii<4; ++ii) {
  474. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii); //[ipc] ;
  475. C(ii, 0) = 1;
  476. C(ii, 1) = pts[0];
  477. C(ii, 2) = pts[1];
  478. C(ii, 3) = pts[2];
  479. }
  480. Real V = (1./6.) * C.determinant(); // volume of tetrahedra
  481. //Real V = C.determinant(); // volume of tetrahedra
  482. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  483. int ID[4];
  484. ID[0] = Ids->GetId(0);
  485. ID[1] = Ids->GetId(1);
  486. ID[2] = Ids->GetId(2);
  487. ID[3] = Ids->GetId(3);
  488. /* bad news bears for magnet */
  489. double* pt = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(0);
  490. Real sum(0);
  491. /*
  492. if (!pe) {
  493. if (std::abs(pt[0]) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  494. sum = 1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  495. pe = true;
  496. }
  497. }*/
  498. if (ID[0] == 26) {
  499. sum = 10;
  500. }
  501. if (ID[0] == 30) {
  502. sum = -10;
  503. }
  504. /*
  505. if (!ne) {
  506. if (std::abs(pt[0]+1.) < .2 && std::abs(pt[1]-15) < .2 && pt[2] < 3.25 ) {
  507. sum = -1; //vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0];
  508. std::cout << "Negative Electroce\n";
  509. ne = true;
  510. }
  511. }
  512. */
  513. //for (int ii=0; ii<4; ++ii) {
  514. //g(ID[ii]) += (V/4.) * ( vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0] ) ;
  515. //if ( std::abs(vtkGrid->GetPointData()->GetScalars()->GetTuple(ID[ii])[0]) > 1e-3 )
  516. //}
  517. // TODO check Load Vector...
  518. g(ID[0]) = sum; //(V/4.) * sum; // Why 3, Why!!!?
  519. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfCells()))) << std::flush ;
  520. }
  521. /*
  522. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  523. vtkGrid->GetPoint(ic, r0);
  524. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  525. double g0 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ;
  526. //std::cout << "num conn " << connectedVertices->GetNumberOfIds() << std::endl;
  527. if ( std::abs(g0) > 1e-3 ) {
  528. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  529. int ii = connectedVertices->GetId(i);
  530. vtkGrid->GetPoint(ii, r1);
  531. double g1 = vtkGrid->GetPointData()->GetScalars()->GetTuple(ii)[0] ;
  532. //g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  533. if ( std::abs(g1) > 1e-3 ) {
  534. g(ic) += CompositeSimpsons2(g1, g0, r1, r0, 1000);
  535. }
  536. //g(ic) += CompositeSimpsons2(g0, r1, r0, 8);
  537. //if ( std::abs(g1) > 1e-3 ) {
  538. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 8);
  539. //g(ic) += CompositeSimpsons2(g0, g1, r0, r1, 100); // / (2*dist(r0,r1)) ;
  540. // g(ic) += g0*dist(r0,r1); //CompositeSimpsons2(g0, r0, r1, 8);
  541. //g(ic) += CompositeSimpsons2(g0, r0, r1, 8);
  542. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  543. //} //else {
  544. // g(ic) += g0; //CompositeSimpsons2(g0, r0, r1, 8);
  545. //}
  546. }
  547. }
  548. //g(ic) = 2.* vtkGrid->GetPointData()->GetScalars()->GetTuple(ic)[0] ; // Why 2?
  549. //std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  550. }
  551. */
  552. } else if (vtkG) { // VTK implicit function, proceed with care
  553. std::cout << "\nUsing implicit file from file" << std::endl;
  554. // OpenMP right here
  555. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  556. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  557. //vtkGrid->GetPoint(ic, r0);
  558. //g(ic) += vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  559. // std::cout << vtkG->FunctionValue(r0[0], r0[1], r0[2]) << std::endl;
  560. //g(ic) += vtkGrid->GetPointData()->GetScalars()->GetTuple1(ic);// FunctionValue(r0[0], r0[1], r0[2]) ;
  561. /*
  562. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  563. int ii = connectedVertices->GetId(i);
  564. vtkGrid->GetPoint(ii, r1);
  565. g(ic) += CompositeSimpsons2(vtkG, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  566. }
  567. */
  568. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  569. }
  570. } else if (gFcn3) {
  571. std::cout << "\nUsing gFcn3" << std::endl;
  572. for (int ic=0; ic < vtkGrid->GetNumberOfPoints(); ++ic) {
  573. vtkSmartPointer<vtkIdList> connectedVertices = GetConnectedPoints(ic);
  574. vtkGrid->GetPoint(ic, r0);
  575. // TODO, test OMP sum reduction here. Is vtkGrid->GetPoint thread safe?
  576. //Real sum(0.);
  577. //#ifdef LEMMAUSEOMP
  578. //#pragma omp parallel for reduction(+:sum)
  579. //#endif
  580. for(vtkIdType i = 0; i < connectedVertices->GetNumberOfIds(); ++i) {
  581. int ii = connectedVertices->GetId(i);
  582. vtkGrid->GetPoint(ii, r1);
  583. g(ic) += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  584. //sum += CompositeSimpsons2(gFcn3, r0, r1, 8); // vtkG->FunctionValue(r0[0], r0[1], r0[2]) ;
  585. }
  586. std::cout << "\r" << (int)(1e2*((float)(ic) / (float)(vtkGrid->GetNumberOfPoints()))) << std::flush ;
  587. //g(ic) = sum;
  588. }
  589. } else {
  590. std::cout << "No source specified\n";
  591. exit(EXIT_FAILURE);
  592. }
  593. // std::cout << g << std::endl;
  594. //g(85) = 1;
  595. std::cout << "\nSolving" << std::endl;
  596. ////////////////////////////////////////////////////////////
  597. // Solving:
  598. //Eigen::SimplicialCholesky<Eigen::SparseMatrix<Real>, Eigen::Lower > chol(A); // performs a Cholesky factorization of A
  599. //VectorXr u = chol.solve(g);
  600. //Eigen::SparseLU<Eigen::SparseMatrix<Real, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > solver;
  601. //solver.analyzePattern(A);
  602. //solver.factorize(A);
  603. //VectorXr u = solver.solve(g);
  604. //Eigen::ConjugateGradient<Eigen::SparseMatrix<Real, Eigen::Lower > Eigen::DiagonalPreconditioner > cg;
  605. Eigen::ConjugateGradient< Eigen::SparseMatrix<Real> > cg(A);
  606. //Eigen::BiCGSTAB<Eigen::SparseMatrix<Real> > cg(A);
  607. cg.setMaxIterations(3000);
  608. //cg.compute(A);
  609. //std::cout << "Computed " << std::endl;
  610. VectorXr u = cg.solve(g);
  611. std::cout << "#iterations: " << cg.iterations() << std::endl;
  612. std::cout << "estimated error: " << cg.error() << std::endl;
  613. vtkDoubleArray *gArray = vtkDoubleArray::New();
  614. vtkDoubleArray *uArray = vtkDoubleArray::New();
  615. uArray->SetNumberOfComponents(1);
  616. gArray->SetNumberOfComponents(1);
  617. for (int iu = 0; iu<u.size(); ++iu) {
  618. uArray->InsertTuple1(iu, u[iu]);
  619. gArray->InsertTuple1(iu, g[iu]);
  620. }
  621. uArray->SetName("u");
  622. gArray->SetName("g");
  623. vtkGrid->GetPointData()->AddArray(uArray);
  624. vtkGrid->GetPointData()->AddArray(gArray);
  625. vtkXMLDataSetWriter *Writer = vtkXMLDataSetWriter::New();
  626. Writer->SetInputData(vtkGrid);
  627. Writer->SetFileName(fname.c_str());
  628. Writer->Write();
  629. Writer->Delete();
  630. Surface->Delete();
  631. gArray->Delete();
  632. uArray->Delete();
  633. }
  634. // Uses simpon's rule to perform a definite integral of a
  635. // function (passed as a pointer). The integration occurs from
  636. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  637. Real FEM4EllipticPDE::CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int n) {
  638. Vector3r R0(r0[0], r0[1], r0[2]);
  639. Vector3r R1(r1[0], r1[1], r1[2]);
  640. // force n to be even
  641. assert(n > 0);
  642. //n += (n % 2);
  643. Real h = dist(r0, r1) / (Real)(n) ;
  644. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  645. Vector3r dr = (R1 - R0).array() / Real(n);
  646. Vector3r rx;
  647. rx.array() = R0.array() + dr.array();
  648. for (int i=1; i<n; i+=2) {
  649. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]);
  650. rx += 2.*dr;
  651. }
  652. rx.array() = R0.array() + 2*dr.array();
  653. for (int i=2; i<n-1; i+=2) {
  654. S += 2.*f->FunctionValue(rx[0], rx[1], rx[2]);
  655. rx += 2.*dr;
  656. }
  657. return h * S / 3.;
  658. }
  659. // Uses simpon's rule to perform a definite integral of a
  660. // function (passed as a pointer). The integration occurs from
  661. // (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  662. // This is just here as a convenience
  663. Real FEM4EllipticPDE::CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int n) {
  664. return dist(r0,r1)*f;
  665. /*
  666. Vector3r R0(r0[0], r0[1], r0[2]);
  667. Vector3r R1(r1[0], r1[1], r1[2]);
  668. // force n to be even
  669. assert(n > 0);
  670. //n += (n % 2);
  671. Real h = dist(r0, r1) / (Real)(n) ;
  672. Real S = f + f;
  673. Vector3r dr = (R1 - R0).array() / Real(n);
  674. //Vector3r rx;
  675. //rx.array() = R0.array() + dr.array();
  676. for (int i=1; i<n; i+=2) {
  677. S += 4. * f;
  678. //rx += 2.*dr;
  679. }
  680. //rx.array() = R0.array() + 2*dr.array();
  681. for (int i=2; i<n-1; i+=2) {
  682. S += 2. * f;
  683. //rx += 2.*dr;
  684. }
  685. return h * S / 3.;
  686. */
  687. }
  688. /*
  689. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  690. * test function owned by the FEM implimentaion.
  691. */
  692. Real FEM4EllipticPDE::CompositeSimpsons2(vtkImplicitFunction* f,
  693. Real r0[3], Real r1[3], int n) {
  694. Vector3r R0(r0[0], r0[1], r0[2]);
  695. Vector3r R1(r1[0], r1[1], r1[2]);
  696. // force n to be even
  697. assert(n > 0);
  698. //n += (n % 2);
  699. Real h = dist(r0, r1) / (Real)(n) ;
  700. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  701. Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  702. //Real S = f->FunctionValue(r0) + f->FunctionValue(r1);
  703. Vector3r dr = (R1 - R0).array() / Real(n);
  704. Vector3r rx;
  705. rx.array() = R0.array() + dr.array();
  706. for (int i=1; i<n; i+=2) {
  707. S += 4. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  708. rx += 2.*dr;
  709. }
  710. rx.array() = R0.array() + 2*dr.array();
  711. for (int i=2; i<n-1; i+=2) {
  712. S += 2. * f->FunctionValue(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  713. rx += 2.*dr;
  714. }
  715. return h * S / 3.;
  716. }
  717. /*
  718. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  719. * test function owned by the FEM implimentaion.
  720. */
  721. Real FEM4EllipticPDE::CompositeSimpsons2( Real (*f)(const Real&, const Real&, const Real&),
  722. Real r0[3], Real r1[3], int n) {
  723. Vector3r R0(r0[0], r0[1], r0[2]);
  724. Vector3r R1(r1[0], r1[1], r1[2]);
  725. // force n to be even
  726. assert(n > 0);
  727. //n += (n % 2);
  728. Real h = dist(r0, r1) / (Real)(n) ;
  729. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  730. //Real S = f(r0[0], r0[1], r0[2])*Hat(R0, R0, R1) + f(r1[0], r1[1], r1[2])*Hat(R1, R0, R1);
  731. Real S = f(r0[0], r0[1], r0[2]) + f(r1[0], r1[1], r1[2]);
  732. Vector3r dr = (R1 - R0).array() / Real(n);
  733. Vector3r rx;
  734. rx.array() = R0.array() + dr.array();
  735. for (int i=1; i<n; i+=2) {
  736. S += 4. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  737. rx += 2.*dr;
  738. }
  739. rx.array() = R0.array() + 2*dr.array();
  740. for (int i=2; i<n-1; i+=2) {
  741. S += 2. * f(rx[0], rx[1], rx[2]) * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  742. rx += 2.*dr;
  743. }
  744. return h * S / 3.;
  745. }
  746. /*
  747. * Performs numerical integration of two functions, one is constant valued f, the other is the FEM
  748. * test function owned by the FEM implimentaion.
  749. */
  750. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int n) {
  751. Vector3r R0(r0[0], r0[1], r0[2]);
  752. Vector3r R1(r1[0], r1[1], r1[2]);
  753. // force n to be even
  754. assert(n > 0);
  755. //n += (n % 2);
  756. Real h = dist(r0, r1) / (Real)(n) ;
  757. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  758. Real S = 2*f; //*Hat(R0, R0, R1) + f*Hat(R1, R0, R1);
  759. Vector3r dr = (R1 - R0).array() / Real(n);
  760. Vector3r rx;
  761. rx.array() = R0.array() + dr.array();
  762. for (int i=1; i<n; i+=2) {
  763. S += 4. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  764. rx += 2.*dr;
  765. }
  766. rx.array() = R0.array() + 2*dr.array();
  767. for (int i=2; i<n-1; i+=2) {
  768. S += 2. * f * Hat(rx, R0, R1) * Hat(rx, R1, R0);
  769. rx += 2.*dr;
  770. }
  771. return h * S / 3.;
  772. }
  773. /*
  774. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  775. * test function owned by the FEM implimentaion.
  776. */
  777. Real FEM4EllipticPDE::CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  778. Vector3r R0(r0[0], r0[1], r0[2]);
  779. Vector3r R1(r1[0], r1[1], r1[2]);
  780. // force n to be even
  781. assert(n > 0);
  782. //n += (n % 2);
  783. Real h = dist(r0, r1) / (Real)(n) ;
  784. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  785. // NO! We are looking at 1/2 hat often!!! So S = f1!
  786. Real S = f1; //f0*Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  787. Vector3r dr = (R1 - R0).array() / Real(n);
  788. // linear interpolate function
  789. //Vector3r rx;
  790. //rx.array() = R0.array() + dr.array();
  791. for (int i=1; i<n; i+=2) {
  792. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  793. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0) ;
  794. }
  795. //rx.array() = R0.array() + 2*dr.array();
  796. for (int i=2; i<n-1; i+=2) {
  797. double fx = f0 + (f1 - f0) * ((i*h)/(h*n));
  798. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1);// * Hat(R1.array() + i*dr.array(), R1, R0);
  799. }
  800. return h * S / 3.;
  801. }
  802. /*
  803. * Performs numerical integration of two functions, one is passed as a vtkImplicitFunction, the other is the FEM
  804. * test function owned by the FEM implimentaion.
  805. */
  806. Real FEM4EllipticPDE::CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int n) {
  807. Vector3r R0(r0[0], r0[1], r0[2]);
  808. Vector3r R1(r1[0], r1[1], r1[2]);
  809. // force n to be even
  810. assert(n > 0);
  811. //n += (n % 2);
  812. Real h = dist(r0, r1) / (Real)(n) ;
  813. // For Gelerkin (most) FEM, we can omit this one as test functions are zero valued at element boundaries
  814. // NO! We are looking at 1/2 hat often!!! So S = f1!
  815. Real S = f0+f1; //Hat(R0, R0, R1) + f1*Hat(R1, R0, R1);
  816. Vector3r dr = (R1 - R0).array() / Real(n);
  817. // linear interpolate function
  818. //Vector3r rx;
  819. //rx.array() = R0.array() + dr.array();
  820. for (int i=1; i<n; i+=2) {
  821. double fx = 1;// f0 + (f1 - f0) * ((i*h)/(h*n));
  822. S += 4. * fx * Hat(R0.array() + i*dr.array(), R0, R1) * Hat(R1.array() + i*dr.array(), R1, R0) ;
  823. }
  824. //rx.array() = R0.array() + 2*dr.array();
  825. for (int i=2; i<n-1; i+=2) {
  826. double fx = 1; // f0 + (f1 - f0) * ((i*h)/(h*n));
  827. S += 2. * fx * Hat(R0.array() + i*dr.array(), R0, R1)* Hat(R1.array() + i*dr.array(), R1, R0);
  828. }
  829. return h * S / 3.;
  830. }
  831. //--------------------------------------------------------------------------------------
  832. // Class: FEM4EllipticPDE
  833. // Method: Hat
  834. //--------------------------------------------------------------------------------------
  835. Real FEM4EllipticPDE::Hat ( const Vector3r& r, const Vector3r& r0, const Vector3r& r1 ) {
  836. //return (r-r0).norm() / (r1-r0).norm() ;
  837. return dist(r, r0) / dist(r1, r0) ;
  838. } // ----- end of method FEM4EllipticPDE::Hat -----
  839. } // ----- end of Lemma name -----