Galerkin FEM for elliptic PDEs
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

FEM4EllipticPDE.cpp 46KB

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