Hi,

thank you for your input. I will think a bit... If I come up with something 
useful I'll let you know




Joa

On Tue, May 30, 2017 at 11:14:10AM -0400, Paul T. Bauman wrote:
> Just FYI, this is something Roy and I have funding to work on over the next
> 2.5 years. We will keep the community posted as this capability comes
> online.
> 
> On Tue, May 30, 2017 at 10:59 AM, John Peterson <jwpeter...@gmail.com>
> wrote:
> 
> > On Tue, May 30, 2017 at 3:29 AM, Joa Ljungvall <libm...@joa.me.uk> wrote:
> >
> > > Dear all,
> > >
> > > I'm revisting an old problem involving solving laplace and poision for a
> > > moderatly complicated geometry. Here I would like to start with a 3D mesh
> > > with
> > > few nodes and refine based on a kerry estimator. I've tried this before
> > > (2009-2010) with libmesh and ran into the problem that I could not figure
> > > out a good way to snap new nodes to the boundary of my geometry. Last
> > time
> > > around I gave up since I created meshes that had all types of problems.
> > >
> > > Now my question, If I make my own MeshRefinment class and rewrite the
> > >
> > > Node * libMesh::MeshRefinement::add_node ( Elem &  parent,
> > >                 unsigned int  child,
> > >                 unsigned int  node,
> > >                 processor_id_type  proc_id
> > >         )
> > > {
> > >
> > > ....
> > > }
> > >
> > > method, more in particular the lines
> > >
> > >   135   // Otherwise we need to add a new node, with a default id and the
> > >   136   // requested processor_id.  Figure out where to add the point:
> > >   137
> > >   138   Point p; // defaults to 0,0,0
> > >   139
> > >   140   for (unsigned int n=0; n != parent.n_nodes(); ++n)
> > >   141     {
> > >   142       // The value from the embedding matrix
> > >   143       const float em_val = parent.embedding_matrix(child,node,n);
> > >   144
> > >   145       if (em_val != 0.)
> > >   146         {
> > >   147           p.add_scaled (parent.point(n), em_val);
> > >   148
> > >   149           // If we'd already found the node we shouldn't be here
> > >   150           libmesh_assert_not_equal_to (em_val, 1);
> > >   151         }
> > >   152     }
> > >
> > > so that if the created node is on a surface I snap it to the boundary
> > > by moving it along the normal of the surface of the parent element. Could
> > > that
> > > work without any strange side effects? This is the old thread of
> > > discussion we
> > > had at the time
> > >
> > > https://sourceforge.net/p/libmesh/mailman/libmesh-users/
> > > thread/20091014130140.GM20082%40joa.me.uk/#msg23751829
> >
> >
> >
> > I don't think deriving your own MeshRefinement object is the right
> > approach... not saying that it couldn't be made to work, just that the
> > class was not really written with that use case in mind, so you would
> > probably run into a lot of gotchas.
> >
> > If it's possible for you to use quadrilateral/hexahedral elements, I think
> > you might have better luck doing what you are trying to do using deal.II.
> >
> > Since the 2009 email thread, they have added support for linking deal.II
> > with OpenCascade [0], including IIRC the ability to snap adaptively refined
> > elements' nodes to the underlying "true" geometry.
> >
> > [0]: https://www.dealii.org/8.5.0/external-libs/opencascade.html
> >
> > --
> > John
> > ------------------------------------------------------------
> > ------------------
> > Check out the vibrant tech community on one of the world's most
> > engaging tech sites, Slashdot.org! http://sdm.link/slashdot
> > _______________________________________________
> > Libmesh-users mailing list
> > Libmesh-users@lists.sourceforge.net
> > https://lists.sourceforge.net/lists/listinfo/libmesh-users
> >

------------------------------------------------------------------------------
Check out the vibrant tech community on one of the world's most
engaging tech sites, Slashdot.org! http://sdm.link/slashdot
_______________________________________________
Libmesh-users mailing list
Libmesh-users@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/libmesh-users

Reply via email to