Go to the documentation of this file.
   23 #include "libmesh/dof_map.h" 
   24 #include "libmesh/elem.h" 
   25 #include "libmesh/enum_xdr_mode.h" 
   26 #include "libmesh/error_vector.h" 
   27 #include "libmesh/equation_systems.h" 
   28 #include "libmesh/explicit_system.h" 
   29 #include "libmesh/libmesh_logging.h" 
   30 #include "libmesh/mesh_base.h" 
   31 #include "libmesh/numeric_vector.h" 
   33 #include "libmesh/exodusII_io.h" 
   34 #include "libmesh/gmv_io.h" 
   35 #include "libmesh/nemesis_io.h" 
   36 #include "libmesh/tecplot_io.h" 
   37 #include "libmesh/xdr_io.h" 
   48   LOG_SCOPE (
"minimum()", 
"ErrorVector");
 
   50   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
   56       libmesh_assert_greater_equal ((*
this)[i], 0.);
 
   58         min = std::min (min, (*
this)[i]);
 
   62   libmesh_assert_greater_equal (min, 0.);
 
   71   LOG_SCOPE (
"mean()", 
"ErrorVector");
 
   73   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
   81         the_mean += ( static_cast<Real>((*
this)[i]) - the_mean ) / (nnz + 1);
 
   94   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
  108       sv.push_back((*
this)[i]);
 
  128   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
  130   LOG_SCOPE (
"variance()", 
"ErrorVector");
 
  132   Real the_variance = 0;
 
  138         const Real delta = ( static_cast<Real>((*
this)[i]) - mean_in );
 
  139         the_variance += (delta * delta - the_variance) / (nnz + 1);
 
  152   LOG_SCOPE (
"cut_below()", 
"ErrorVector");
 
  154   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
  156   std::vector<dof_id_type> cut_indices;
 
  157   cut_indices.reserve(n/2);  
 
  162         if ((*
this)[i] < cut)
 
  164             cut_indices.push_back(i);
 
  176   LOG_SCOPE (
"cut_above()", 
"ErrorVector");
 
  178   const dof_id_type n = cast_int<dof_id_type>(this->size());
 
  180   std::vector<dof_id_type> cut_indices;
 
  181   cut_indices.reserve(n/2);  
 
  186         if ((*
this)[i] > cut)
 
  188             cut_indices.push_back(i);
 
  199   libmesh_assert_less (i, this->size());
 
  206     return ((*
this)[i] != 0.);
 
  213   std::unique_ptr<MeshBase> meshptr = oldmesh.
clone();
 
  218   mesh.allow_renumbering(
false);
 
  219   mesh.all_first_order();
 
  221 #ifdef LIBMESH_ENABLE_AMR 
  224   for (
auto & elem : 
mesh.element_ptr_range())
 
  227       elem->set_p_level(0);
 
  229 #endif // LIBMESH_ENABLE_AMR 
  238   std::vector<dof_id_type> dof_indices;
 
  240   for (
const auto & elem : 
mesh.active_local_element_ptr_range())
 
  250       libmesh_assert_less (elem_id, (*this).size());
 
  254       error_system.
solution->set(solution_index, (*
this)[elem_id]);
 
  262   if (
mesh.max_elem_id() != 
mesh.n_elem() ||
 
  263       mesh.max_node_id() != 
mesh.n_nodes())
 
  265       mesh.allow_renumbering(
true);
 
  266       mesh.renumber_nodes_and_elements();
 
  269   if (filename.rfind(
".gmv") < filename.size())
 
  274   else if (filename.rfind(
".plt") < filename.size())
 
  279 #if defined(LIBMESH_HAVE_EXODUS_API) && defined(LIBMESH_HAVE_NEMESIS_API) 
  280   else if ((filename.rfind(
".nem") < filename.size()) ||
 
  281            (filename.rfind(
".n") < filename.size()))
 
  288 #ifdef LIBMESH_HAVE_EXODUS_API 
  289   else if ((filename.rfind(
".exo") < filename.size()) ||
 
  290            (filename.rfind(
".e") < filename.size()))
 
  297   else if (filename.rfind(
".xda") < filename.size())
 
  304   else if (filename.rfind(
".xdr") < filename.size())
 
  314       libMesh::err << 
"Warning: ErrorVector::plot_error currently only" 
  315                    << 
" supports .gmv, .plt, .xdr/.xda, and .exo/.e (if enabled) output;" << std::endl;
 
  316       libMesh::err << 
"Could not recognize filename: " << filename
 
  
void plot_error(const std::string &filename, const MeshBase &mesh) const
Plots a data file, of a type determined by looking at the file extension in filename,...
 
virtual void write(const std::string &base_filename) override
This method implements writing a mesh to a specified file.
 
The Nemesis_IO class implements reading parallel meshes in the Nemesis file format from Sandia Nation...
 
virtual ErrorVectorReal minimum() const override
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
 
MeshIO class used for writing XDR (eXternal Data Representation) and XDA mesh files.
 
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
 
The libMesh namespace provides an interface to certain functionality in the library.
 
void write_element_data(const EquationSystems &es)
Write out element solution.
 
virtual Real median() override
 
virtual void write(const std::string &) override
This method implements writing a mesh to a specified file.
 
MeshBase * _mesh
Pointer to the mesh, which may be used to decide which elements are active.
 
This class implements writing meshes in the GMV format.
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
virtual const Elem * elem_ptr(const dof_id_type i) const =0
 
This is the MeshBase class.
 
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
 
void write_discontinuous_gmv(const std::string &name, const EquationSystems &es, const bool write_partitioning, const std::set< std::string > *system_names=nullptr) const
Writes a GMV file with discontinuous data.
 
virtual void init()
Initialize all the systems.
 
virtual std::vector< dof_id_type > cut_above(Real cut) const override
 
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
 
bool is_active_elem(dof_id_type i) const
Utility function to decide whether element i is active.
 
virtual std::unique_ptr< MeshBase > clone() const =0
Virtual "copy constructor".
 
This is the EquationSystems class.
 
void write_element_data(const EquationSystems &es)
Write out element solution in parallel, without localizing the solution vector.
 
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
 
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
This class handles the numbering of degrees of freedom on a mesh.
 
This class implements writing meshes in the Tecplot format.
 
virtual Real mean() const override
 
virtual std::vector< dof_id_type > cut_below(Real cut) const override
 
void write(const std::string &name, const XdrMODE, const unsigned int write_flags=(WRITE_DATA), bool partition_agnostic=true) const
Write the systems to disk using the XDR data format.
 
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
 
const DofMap & get_dof_map() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
virtual Real variance() const override