From: Andrew E S. <and...@gm...> - 2012-06-06 21:27:32
|
Program works if I use a wrapper function instead of AnalyticFunction, perhaps I am not using it correctly or there is a bug. The code that works is below, the AnalyticFunction stuff that, at least as I understand, should also work is commented out. If anyone has any information on why this problem exists I would appreciated it, I prefer using AnalyticFunction if I can. Peace, Andrew --- PROGRAM --- // Standard library includes #include <iostream> #include <algorithm> #include <math.h> #include <stdio.h> // libMesh includes #include "libmesh.h" #include "mesh.h" #include "mesh_generation.h" #include "equation_systems.h" #include "transient_system.h" #include "explicit_system.h" #include "analytic_function.h" #include "vtk_io.h" #include "exodusII_io.h" // Bring in everything from the libMesh namespace using namespace libMesh; // A function for initilizing thermal conductivity #include <boost/math/constants/constants.hpp> Number initial_conductivity (const Point& p, const Real t = 0){ const double pi = boost::math::constants::pi<double>(); return 1 + exp(-t) * sin(pi*p(0)) * sin(pi*p(1)); } // Wrapper function for exact solution, used for initilization Number initial_function(const Point& p, const Parameters&, const std::string&, const std::string&){ return initial_conductivity(p, 0.0); } // Begin main function int main (int argc, char** argv){ // Initialize libraries LibMeshInit init (argc, argv); // Generate a mesh Mesh mesh; MeshTools::Generation::build_square(mesh, 10, 10, 0., 1, 0., 1, QUAD8); mesh.all_second_order(); // Create a system for storing nodal values of conductivity EquationSystems eq_sys(mesh); eq_sys.add_system<TransientExplicitSystem>("data"); eq_sys.get_system("data").add_variable("k", SECOND, MONOMIAL); eq_sys.init(); // Project the initial values of conductivity //eq_sys.get_system("data").time = 0; //AnalyticFunction<Number> func_object(initial_conductivity); //eq_sys.get_system("data").project_solution(&func_object); eq_sys.get_system("data").project_solution(initial_function, NULL, eq_sys.parameters); // Write conductivity values to a file ExodusII_IO(mesh).write_equation_systems("example3.ex2", eq_sys); return 0; } On Wed, Jun 6, 2012 at 2:28 PM, Andrew E Slaughter < and...@gm...> wrote: > I searched through the mailings lists and found some discussion that > explains that nodal and element data should be stored in an ExplicitSystem, > but could not find any example or explanation to how it should be done. > > So, I am working on a test program, but run into a problem when projecting > a "solution". What I am attempting to do in my test code is to create a > TransientExplicitSystem to store a parameter (conductivity). Then using an > equation, I initialize the value across the mesh via projection (I imagine > that I should be able to use projection at any point in time to alter the > values). > > The program compiles and links, but produces a PETSc segmentation fault > when executed. I don't understand the concept of using systems for data > storing so I don't have a good feel for what I might be doing wrong. > > I would greatly appreciate any help getting this up and running or any > recommendations in general for storing and updating nodal and/or elemental > data. > > Thanks, > Andrew > > --- MY PROGRAM --------- > > // Standard library includes > #include <iostream> > #include <algorithm> > #include <math.h> > #include <stdio.h> > > // libMesh includes > #include "libmesh.h" > #include "mesh.h" > #include "mesh_generation.h" > #include "equation_systems.h" > #include "transient_system.h" > #include "explicit_system.h" > #include "analytic_function.h" > #include "vtk_io.h" > > // Bring in everything from the libMesh namespace > using namespace libMesh; > > // A function for initilizing thermal conductivity > #include <math.h> > #include <boost/math/constants/constants.hpp> > Number initial_conductivity (const Point& p, const Real t){ > const double pi = boost::math::constants::pi<double>(); > return 1 + exp(-t) * sin(pi*p(0)) * sin(pi*p(1)); > } > > // Begin main function > int main (int argc, char** argv){ > > // Initialize libraries > LibMeshInit init (argc, argv); > > // Generate a mesh > Mesh mesh; > MeshTools::Generation::build_square(mesh, 10, 10, 0., 1, 0., 1, QUAD4); > mesh.all_second_order(); > > // Create a system for storing nodal values of conductivity > EquationSystems eq_sys(mesh); > eq_sys.add_system<TransientExplicitSystem>("data"); > eq_sys.get_system("data").add_variable("k", SECOND, MONOMIAL); > eq_sys.init(); > // Project the initial values of conductivity > AnalyticFunction<Number> func_object(initial_conductivity); > eq_sys.get_system("data").project_solution(&func_object); // > segmentation fault here > > // Write conductivity values to a file > //... > > return 0; > } > > -- > Andrew E. Slaughter, PhD > and...@gm... > > Materials Process Design and Control Laboratory > Sibley School of Mechanical and Aerospace Engineering > 169 Frank H. T. Rhodes Hall > Cornell University > Ithaca, NY 14853-3801 > (607) 229-1829 > -- Andrew E. Slaughter, PhD and...@gm... Materials Process Design and Control Laboratory Sibley School of Mechanical and Aerospace Engineering 169 Frank H. T. Rhodes Hall Cornell University Ithaca, NY 14853-3801 (607) 229-1829 |