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

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