Hi again,
during the day I been working a bit on the problem: Currently I have the
exact geometry inside my program (implemented via CGS and GEANT4). What I
have tried is to, after the refinement move the new boundary points using
a code like this one (it tries to move all nodes on a boundary to the correct
"radius" along the rhat vector
void libMESH::AdjustMesh(Mesh &mesh,BaseBoundary *BC)
{
MeshBase::element_iterator el =
mesh.active_elements_begin();
const MeshBase::element_iterator end_el =
mesh.active_elements_end();
for ( ; el != end_el; ++el)
{
for (unsigned int s=0; s<(*el)->n_sides(); s++)
if ((*el)->neighbor(s) == NULL) {
AutoPtr<DofObject> face = (*el)->side(s);
for(int atnode=0; atnode<3; ++atnode){
Node *node = (dynamic_cast<Elem*>(face.get())->set_node(atnode));
//if z==0 its ok
if(fabs((*node)(2))>1e-6){
if((*node)(2)<BC->CylLength()){//Simple case
//Check which boundary is the closest
double Router = BC->RAtZ(MakeHepVector((*node)));
double Rinner =
BC->RadiusOnCentralContact(MakeHepVector((*node)));
double Rnode =
sqrt((*node)(0)*(*node)(0)+(*node)(1)*(*node)(1));
if(fabs(Router-Rnode)<fabs(Rinner-Rnode)){
(*node)(0)*=Router/Rnode;
(*node)(1)*=Router/Rnode;
} else {
(*node)(0)*=Rinner/Rnode;
(*node)(1)*=Rinner/Rnode;
}
} else {//Difficult, we have to check if the node should be on
//one of the "circles"
}
}
}
}
}
}
where BaseBoundary is a class to take care of the geometry. Trying with this I
get
EquationSystems
n_systems()=1
System "Laplace"
Type "LinearImplicit"
Variables="u"
Finite Element Types="LAGRANGE"
Approximation Orders="FIRST"
n_dofs()=715
n_local_dofs()=715
n_constrained_dofs()=0
n_vectors()=1
Beginning Solve 0
System has: 715 degrees of freedom.
Linear solver converged at step: 5000, final residual: 1.1091e-16
Refining the mesh...
Beginning Solve 1
ERROR: negative Jacobian: -83.5253 in element 7069
[0] src/fe/fe_map.C, line 341, compiled Oct 11 2009 at 15:03:39
terminate called after throwing an instance of 'libMesh::LogicError'
what(): Error in libMesh internal logic
Abort trap
It is quite possible that I move the wrong nodes the wrong way, but is my
problem ;). The general question is, is this a possible route to take? I call
this routine after refinement and coursening, before equation_systems.reinit().
cheers
Joa
On Tue, Oct 13, 2009 at 10:54:57AM -0500, John Peterson wrote:
> On Tue, Oct 13, 2009 at 10:36 AM, Joa Ljungvall <[email protected]> wrote:
> > Hi,
> >
> > The problem is that I have to solve my equation on the same domain with
> > very different boundary conditions, making a uniform refinement very
> > inefficient and a waste of time, and worse in my case, RAM. So I do need
> > to start with a very course mesh and refine.
>
> As David explained, there is unfortunately currently no "feedback"
> mechanism between the mesh refinement procedure and an underlying
> geometric description.
>
> One approach might be to maintain a suitably fine, lower-dimensional
> (and hence, low memory) discretization of the boundary as a secondary
> Mesh object in LibMesh, and then somehow use this fine discretization
> when placing new points, but again, there is no existing mechanism for
> doing this currently in the library.
>
> --
> John
------------------------------------------------------------------------------
Come build with us! The BlackBerry(R) Developer Conference in SF, CA
is the only developer event you need to attend this year. Jumpstart your
developing skills, take BlackBerry mobile applications to market and stay
ahead of the curve. Join us from November 9 - 12, 2009. Register now!
http://p.sf.net/sfu/devconference
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users