>> Patch didn't apply cleanly for me (submit as an attachment next time?
>> my email system must have mangled the whitespace or something) but it
>> looks simple enough to apply by hand.
Sorry about that. I'll attach patches as files in the future.
>
> And I'm glad I looked at it more closely - shouldn't we also be using
> vtkDoubleArray in system_vectors_to_vtk? It looks as if your patch
> bumps up the precision on system.solution but not on the other vectors
> being output.
Yes, that would seem like a good idea to me. However,
system_vectors_to_vtk is not actively being used (the only call is
commented and even bracketed by #ifdef DEBUG), and frankly, it is a
huge mess. I suggest commenting that method out, with a remark that it
doesn't support AMR or higher order variables (it really doesn't -- it
relies on n_dof() == n_nodes(), which is fatal). Besides, I don't
think that VTKIO::write_equation_systems should write system vectors
anyway. Neither of ExodusII_IO or GMVIO do.
> Does VTK have native support for higher-order data? If so, adding
> that to our VTKIO output would be fantastic.
Yes, VTK does, but it's not implemented in VTKIO. I attached a patch
that does it. :-)
Note: The patch only modifies solution_to_vtk, while leaving
system_vectors_to_vtk untouched. See above for the reason why...
> If not, then it ought to be possible to implement a more general
> solution. "Subdivide an order-p element into p^d order-1 elements" is
> a pretty file-format-independent task, even if we limited the initial
> testing to VTKIO.
Agreed. It's just that VTK output is implemented fundamentally
different than others. But it should be relatively easy to modify
EquationSystems::build_solution_vector and
EquationSystems::build_variable_names, which seem to be only used by
ExodusII_IO and GMVIO currently, in the same way, i.e. such that
build_variable_names also returns all components with a suffix, and
build_solution_vector returns also the higher order components. It
seems to me that this would be sufficient already.
Roman
Index: src/mesh/vtk_io.C
===================================================================
--- src/mesh/vtk_io.C (revision 4200)
+++ src/mesh/vtk_io.C (working copy)
@@ -303,9 +303,10 @@
// write only single processor
if(libMesh::processor_id()==0){
- const MeshBase& mesh = (MeshBase&)es.get_mesh();
+ const MeshBase& mesh = (MeshBase&)es.get_mesh();
+ const unsigned int n_nodes = mesh.n_nodes();
+ libmesh_assert(n_nodes > 0);
- const unsigned int n_nodes = es.get_mesh().n_nodes();
// loop over the systems
for(unsigned int i=0;i<es.n_systems();++i){ // for all systems, regardless of whether they are active or not
@@ -315,24 +316,25 @@
// loop over variables
for(unsigned int j=0;j<n_vars;++j){
- std::string name = sys.variable_name(j);
+ const unsigned int n_comp = mesh.node(0).n_comp(i,j);
+ std::string name = sys.variable_name(j);
+ vtkDoubleArray *data = vtkDoubleArray::New();
+ data->SetName(name.c_str());
+ data->SetNumberOfComponents(n_comp);
+ data->SetNumberOfTuples(n_nodes);
- vtkDoubleArray *data = vtkDoubleArray::New();
+ for(unsigned int k=0;k<n_nodes;++k){
- data->SetName(name.c_str());
+ for(unsigned int c=0;c<n_comp;++c){
- data->SetNumberOfValues(n_nodes);
+ const unsigned int dof_nr = mesh.node(k).dof_number(i,j,c);
- for(unsigned int k=0;k<n_nodes;++k){
-
- const unsigned int dof_nr = mesh.node(k).dof_number(i,j,0);
-
#ifdef LIBMESH_USE_COMPLEX_NUMBERS
- data->SetValue(k,sys.current_solution(dof_nr).real());
+ data->SetComponent(k,c,sys.current_solution(dof_nr).real());
#else
- data->SetValue(k,sys.current_solution(dof_nr));
+ data->SetComponent(k,c,sys.current_solution(dof_nr));
#endif
-
+ }
}
grid->GetPointData()->AddArray(data);
data->Delete();
------------------------------------------------------------------------------
Index, Search & Analyze Logs and other IT data in Real-Time with Splunk
Collect, index and harness all the fast moving IT data generated by your
applications, servers and devices whether physical, virtual or in the cloud.
Deliver compliance at lower cost and gain new business insights.
Free Software Download: http://p.sf.net/sfu/splunk-dev2dev
_______________________________________________
Libmesh-devel mailing list
Libmesh-devel@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/libmesh-devel