From: Roy S. <roy...@ic...> - 2017-11-27 21:23:00
|
So sorry nobody got back to you on this already! On Tue, 14 Nov 2017, Michael Povolotskyi wrote: > I would like to output a nodal data stored in NumericVector, using > > void libMesh::MeshOutput< MT >::write_nodal_data ( const std::string & > fname, > const NumericVector< Number > & parallel_soln, > const std::vector< std::string > & names > ) > > I need to fill NumericVector with values. Well, it needs to be filled with values somehow. Typically those values come from an interpolation/projection of some function or from a solution of some system of equations, rather than being set by hand. > As far as I understand, I need to provide data both for active and non-active > nodes, because the NumericVector contains entities for all local nodes. Oh, it's much more complicated than that: In multiphysics problems, you must have set values for each degree of freedom (DoF) of each variable, and there will therefore in general be multiple variable DoFs per node. In parallel, each processor must have set values for its own local DoFs. On non-Lagrange elements (e.g. hierarchic or discontinuous) or on Lagrange non-isoparametric elements (e.g. a bilinear function on a QUAD9) the common belief that one node == one DoF is false. You may have nodes where a variable has no DoFs (e.g. the mid-edge nodes in that bilinear-on-QUAD9 case), you may have nodes where a variable has multiple DoFs (e.g. the mid-edge of nodes in the hierarchic case for cubic and higher polynomial degree), and you may have DoFs stored on elements where there is no node (e.g. for "bubble function" coefficients for cubic and higher polynomial degree). If you have a function you're trying to discretize, I suggest using System::project_vector() to do the discretization for you. If you insist on doing the interpolation yourself (which isn't completely insane, if you know you're doing Lagrange FE only and you want to be as fast as possible), then the method you want to use to look up DoF indices is DofObject::dof_number(). Give it a system number, variable number, and component number (the latter is always 0 for Lagrange FE on nodes where those FE have a DoF) and it will give you the corresponding DoF index. > 1. Will I get correct order if I iterate from mesh.local_nodes_begin() to > mesh.local_nodes_end() ? Maybe, if the phase of the moon is correct? Even if the answer happens to be yes, definitely don't count on it; it's not a behavior we guarantee. We do guarantee that this iteration will give you the correct *set* of nodes, though, in parallel - local DoFs are stored on and only on local nodes. > 2. Will the values I set for non active nodes matter for the output? Apparently we used to have such a thing as non-active nodes? Today I Learned. These days we strip nodes from the mesh as soon as they're no longer used by any element (if you actually *want* a disconnected node then you create a "NodeElem" to connect it to), so you should never encounter one. --- Roy |