Hi Tom,
Its great that you're getting somewhere with this!
I have taken some good advice from Oded (it turns out we work one building
> over from each other!),
>
What a happy coincidence :-)
> I'm trying hard to maintain the MappingQEulerian approach in my current
> model setup, and I wanted to return to the outline you put forward to see
> if 'Mesquite' can simply be replaced with 'Cubit' in your writeup:
>
In theory this shouldn't be a problem. This is sort of a "black-box" step.
> (exporting from cubit is actually a bit of a problem, I have submitted a
> service email to them - but if I do this through the gui it all works)
>
Here's a snippet of one of my cubit scripts that may be of use?
# === Export geometry ===
# - Cubit format
cubit.cmd('save as "' + out_mesh.output_cub_file + '" overwrite')
# - Abaqus format
cubit.cmd('set Abaqus precision 6')
cubit.cmd('export Abaqus "' + out_mesh.output_abaqus_file + '" dimension '
+ str(in_geom.dim) + ' overwrite')
2. A function that actually
>> <https://www.dealii.org/developer/doxygen/deal.II/step_18.html#TopLevelmove_mesh>
>>
>> moves
>> <https://www.dealii.org/developer/doxygen/deal.II/step_42.html#PlasticityContactProblemmove_mesh>
>>
>> the triangulation vertices for you based on this vector.
>>
>
> At this point the Triangulation attached to my GridIn object is the new
> triangulation, so this may not be a problem.
>
Ok, so this might be a point of difficulty. See the answer to your next
question...
>
>
>> 3. A function that would interpolate the optimal Q1 Euler vector onto a
>> given FE space, which would presumably represent the displacement solution
>> to some other problem.
>>
>>
> *Here's where my question for you comes in:* At this stage I have a new
> Triangulation (call it temp_triangulation) which holds my optimal
> vertices. The vertices in temp_triangulation live on some surface in real
> space. I have my original Triangulation (call it orig_triangulation) which
> satisfies orig_triangulation.n_vertices() ==
> temp_triangulation.n_vertices(). *How to I obtain a new euler_vector of
> size dof_handler.n_dofs() at this stage? *
>
> I have done the following so far:
>
> 1 /* update euler_vector to describe smoothed mesh */
> 2 FESystem<dim,spacedim> vertex_fe(FE_Q<dim,spacedim>(1), spacedim);
> 3 DoFHandler<dim,spacedim> vertex_dof_handler(temp_triangulation);
> 4 vertex_dof_handler.distribute_dofs(vertex_fe);
> 5 Vector<double> vertex_euler_vector(vertex_dof_handler.n_dofs());
> 6
> 7 std::vector<Point<spacedim> > original_vertices =
> triangulation.get_vertices();
> 8 std::vector<bool>
> vertex_touched(temp_triangulation.n_vertices(),false);
> 9
> 10 typename DoFHandler<dim,spacedim>::active_cell_iterator
> 11 cell=vertex_dof_handler.begin_active(),
> 12 endc=vertex_dof_handler.end();
> 13 for (; cell!=endc; ++cell)
> 14 for (unsigned int v=0; v < 4; ++v)
>
>
By "4", you mean GeometryInfo<spacedim>::vertices_per_cell, right?
> 15 if (!vertex_touched[cell->vertex_index(v)])
> 16 {
> 17 for (unsigned int k=0; k<spacedim; ++k)
> 18 {
> 19 vertex_euler_vector(cell->vertex_dof_index(v,k))
> 20 = cell->vertex(v)(k) -
> original_vertices[cell->vertex_index(v)](k);
> 21 }
> 22 vertex_touched[cell->vertex_index(v)] = true;
> 23 }
>
>
>
At first glance this looks ok. My biggest concern is that you've assumed
that your old and new triangulations have the same vertex ordering. I'm
inclined to think that this is ok because I don't think that Cubit
renumbers mesh entities unless you ask it to, and nor does deal.II. But you
should probably hand-check this on a small geometry.
Basically, I have built vertex_euler_vector to map the vertices of my
> original reference domain to new points which are more or less on my
> original surface. This holds essentially Q1 information - only information
> at the nodes of the original reference domain. I *think* that what I
> need to do now is interpolate vertex_euler_vector onto my original finite
> element space (which was Q2, but could be anything). If that is what is
> needed, I would appreciate some advice on how to accomplish it!
>
I think that VectorTools::interpolate_to_different_mesh
<https://www.dealii.org/developer/doxygen/deal.II/namespaceVectorTools.html#af68148d58c8dfd0916eceab9d89d74d5>
could
do the job for you; even though it needs the same FE discretisation I
*think* this means that both DoFHandlers (both based on the original
triangulation) must be structured like FESystem<dim>(FE_Q(n),dim). If I'm
wrong on this then you could wrap your Euler vector in an FEFieldFunction
<https://www.dealii.org/developer/doxygen/deal.II/classFunctions_1_1FEFieldFunction.html>
and just use apply of the VectorTools::interpolate
<https://www.dealii.org/developer/doxygen/deal.II/namespaceVectorTools.html#a05db6c8cebf924b417dd92f525efe3db>
functions - this may be slow if the triangulation is huge, but it would be
a good first step to take.
> Does any of this sound reasonable?
>
Well, none of it sounds crazy :-)
J-P
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
For more options, visit https://groups.google.com/d/optout.