|
From: Derek G. <fri...@gm...> - 2006-06-22 19:55:28
|
Ok got it... here's a stack trace on a failing call: ########### #0 dumb_ass () at fe_map.C:34 #1 0x00002afe1d47fd5a in FE<2u, (libMeshEnums::FEFamily)0>::inverse_map (elem=0xed52b0, physical_point=@0xfe4db0, tolerance=0.0001, secure=true) at fe_map.C:1364 #2 0x00002afe1d45b565 in FEInterface::inverse_map (dim=2, fe_t=@0xeec470, elem=0xed52b0, p=@0xfe4db0, tolerance=0.0001, secure=true) at fe_interface.C:1186 #3 0x00002afe1d5f7d6c in MeshFunction::operator() (this=0x7fffffc8b2d0, p=@0xfe4db0, output=@0x7fffffc8b080) at mesh_function.C:190 #4 0x00002afe1d5f7c2c in MeshFunction::operator() (this=0x7fffffc8b2d0, p=@0xfe4db0, time=0) at mesh_function.C:157 #5 0x000000000044a155 in main (argc=5, argv=0x7fffffc8ba68) at ex10.C:461 ############ Note that I had to put an un-templated function into fe_map.C... then I called it just before the WARNING was printed... I was able to put a break point on that function therefore finally recovering the stack at one that fails. The problem was that gdb couldn't figure out which templated version of inverse_map to put a break on... and it was always guess wrong and therefore never stopping... very weird. But I worked around it ;-) The "problem" appears to be happening when mesh_function is looking for the mapped_point. I don't know enough about this stuff to interpret though. Let me know if you guys need me to do anything else. Derek On 6/21/06, Roy Stogner <roy...@ic...> wrote: > On Wed, 21 Jun 2006, Derek Gaston wrote: > > > what I get: > > > > WARNING: Newton scheme has not converged in 31 iterations: > > physical_point=(x,y,z)=(0.158938, 0.422823, 0) > > physical_guess=(x,y,z)=(0.497654, 0.17708, 0) > > dp=(x,y,z)=(-9.77611, 7.05057, 0) > > p=(x,y,z)=(0.685321, 0.746925, 0) > > secure=1 error=12.0533 > > > > How in the _HELL_ is secure=1?!???!???!?!! > > Keep in mind that secure==true is the default argument, and > inverse_map gets called from more than just the point locator. Even > when you're not using adaptivity, your code may be calling inverse_map > from FE::nodal_soln, FE::edge_reinit, or MeshFunction::operator(). > However, none of those functions should ever give you an inverse_map > that doesn't converge, so it's important that we find what's wrong. > > At the very least, it would be a good idea to print out the node > locations of the element for which inverse_map is failing. If the > physical point is inside the element and the element has a positive > lower bound on it's jacobian, then we'll know it's a serious problem > and we'll probably have enough information to replicate it. > --- > Roy > |