25 #include "libmesh/libmesh_base.h"
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/mesh_base.h"
28 #include "libmesh/mesh_output.h"
29 #include "libmesh/equation_systems.h"
30 #include "libmesh/explicit_system.h"
31 #include "libmesh/boundary_mesh.h"
32 #include "libmesh/node.h"
33 #include "libmesh/elem.h"
34 #include "libmesh/enum_order.h"
35 #include "libmesh/enum_elem_type.h"
36 #include "libmesh/numeric_vector.h"
37 #include "libmesh/auto_ptr.h"
38 #include "libmesh/dense_vector.h"
39 #include "libmesh/fe_base.h"
40 #include "libmesh/dof_map.h"
41 #include "libmesh/point.h"
76 void ProjectPointToLine( Point pt, Line line, Real& xi )
78 libmesh_assert( (dim==2) || (dim==3) );
86 for (
unsigned int ii=0; ii<dim; ii++ ) {
87 dx[ii] = p1(ii) - p0(ii);
88 length2 += dx[ii]*dx[ii];
89 xavg[ii] = 0.5*( p1(ii) + p0(ii) );
96 for (
unsigned int ii=0; ii<dim; ii++ ) {
97 xi += 2.0*( pt(ii) - xavg[ii] )*dx[ii];
107 DistanceToSegment( Point pt, Line line, Real& distance )
109 libmesh_assert( (dim==2) || (dim==3) );
112 ProjectPointToLine<dim> (pt, line, xi);
115 if( xi < -1.0 ) xi = -1.0;
116 if( xi > 1.0 ) xi = 1.0;
120 for (
unsigned int ii=0; ii<dim; ii++ ) {
121 xint[ii] = line.p0(ii)*(-0.5*(xi-1.0)) + line.p1(ii)*(0.5*(xi+1.0));
126 for (
unsigned int ii=0; ii<dim; ii++ ) {
127 distance += (pt(ii) - xint[ii])*(pt(ii) - xint[ii]);
129 distance = sqrt(distance);
131 libmesh_assert( (std::isfinite(distance)) && (distance>=0.0) );
149 _equation_systems (es_in),
150 _boundary_mesh (bm_in),
154 libmesh_assert(libMesh::initialized());
160 libMesh::System& sys =
_equation_systems.get_system<libMesh::System>(
"distance_function");
163 sys.add_variable(
"distance", libMesh::FIRST);
166 sys.attach_init_object(*
this);
187 libmesh_assert( node != NULL );
190 libMesh::Real distance = std::numeric_limits<libMesh::Real>::infinity();
194 libmesh_assert( (dim==2) || (dim==3) );
197 std::vector<libMesh::Real> xnode(dim);
198 for(
unsigned int ii=0; ii<dim; ii++ ) xnode[ii] = (*node)(ii);
208 libMesh::MeshBase::const_element_iterator el =
_boundary_mesh.active_elements_begin();
209 const libMesh::MeshBase::const_element_iterator end_el =
_boundary_mesh.active_elements_end();
213 std::cout <<
"There are no elements to iterate through!!! Returning..." << std::endl;
217 const libMesh::Elem* elem_containing_min = NULL;
218 libMesh::Real dmin = std::numeric_limits<Real>::infinity();
220 for ( ; el != end_el; ++el) {
222 const libMesh::Elem* belem = *el;
224 for (
unsigned int bnode=0; bnode<belem->n_nodes(); ++bnode)
226 const libMesh::Point& pbndry = belem->point(bnode);
228 libMesh::Real dnode = 0.0;
229 for (
unsigned int idim=0; idim<dim; ++idim)
230 dnode += (xnode[idim] - pbndry(idim))*(xnode[idim] - pbndry(idim));
232 dnode = std::sqrt(dnode);
237 elem_containing_min = belem;
244 std::set<const libMesh::Elem*> near_min_elems;
245 elem_containing_min->find_point_neighbors (near_min_elems);
248 near_min_elems.insert(elem_containing_min);
250 std::set<const libMesh::Elem*>::iterator nm_el;
253 for (nm_el=near_min_elems.begin(); nm_el!=near_min_elems.end(); ++nm_el)
256 const libMesh::Elem* belem = *nm_el;
259 libmesh_assert( belem->default_order() == FIRST );
264 libMesh::Real dedge = std::numeric_limits<libMesh::Real>::infinity();
274 if (belem->dim()!=1)
continue;
278 libmesh_assert( belem->n_nodes() == 2 );
279 libmesh_assert( belem->type() == EDGE2 );
282 const libMesh::Point& p0 = belem->point(0);
283 const libMesh::Point& p1 = belem->point(1);
289 DistanceToSegment<2>(*node, line, dedge);
291 if ( dedge < distance ) distance = dedge;
302 if (belem->dim()!=2)
continue;
304 libmesh_assert( (belem->type() == TRI3) || (belem->type() == QUAD4) );
306 if ( belem->type() == TRI3 )
314 const libMesh::Point& p0 = belem->point(0);
315 const libMesh::Point& p1 = belem->point(1);
316 const libMesh::Point& p2 = belem->point(2);
318 const libMesh::Real x_xi = (p1(0) - p0(0)), x_et = (p2(0) - p0(0));
319 const libMesh::Real y_xi = (p1(1) - p0(1)), y_et = (p2(1) - p0(1));
320 const libMesh::Real z_xi = (p1(2) - p0(2)), z_et = (p2(2) - p0(2));
322 libMesh::Real A[2][2], b[2];
324 A[0][0] = x_xi*x_xi + y_xi*y_xi + z_xi*z_xi;
325 A[0][1] = x_xi*x_et + y_xi*y_et + z_xi*z_et;
326 A[1][0] = x_xi*x_et + y_xi*y_et + z_xi*z_et;
327 A[1][1] = x_et*x_et + y_et*y_et + z_et*z_et;
329 b[0] = (xnode[0] - p0(0))*x_xi + (xnode[1] - p0(1))*y_xi + (xnode[2] - p0(2))*z_xi;
330 b[1] = (xnode[0] - p0(0))*x_et + (xnode[1] - p0(1))*y_et + (xnode[2] - p0(2))*z_et;
332 const libMesh::Real detA = A[0][0]*A[1][1] - A[1][0]*A[0][1];
333 libmesh_assert( fabs(detA) > 0.0 );
335 const libMesh::Real xi = ( A[1][1]*b[0] - A[0][1]*b[1])/detA;
336 const libMesh::Real et = (-A[1][0]*b[0] + A[0][0]*b[1])/detA;
344 if ( (xi<0.0) || (et<0.0) || (et>1.0-xi) ) {
346 libMesh::Real dtmp = std::numeric_limits<libMesh::Real>::infinity();
354 for (
unsigned int iedge=0; iedge<3; iedge++ ){
359 if (iedge==0) { line.p0 = p1; line.p1 = p2; }
360 else if (iedge==1) { line.p0 = p0; line.p1 = p2; }
361 else { line.p0 = p0; line.p1 = p1; }
363 DistanceToSegment<3>(*node, line, dtmp);
365 if( dtmp < dedge ) dedge = dtmp;
372 libMesh::Real xint[3];
373 for (
unsigned int ii=0; ii<3; ii++ ) {
374 xint[ii] = p0(ii)*(1.0 - xi - et) + p1(ii)*xi + p2(ii)*et;
379 for (
unsigned int ii=0; ii<3; ii++ ) {
380 dedge += (xnode[ii] - xint[ii])*(xnode[ii] - xint[ii]);
387 else if ( belem->type() == QUAD4 )
391 libMesh::Real RTOL = 1e-10;
392 libMesh::Real ATOL = 1e-20;
393 const unsigned int ITER_MAX=1000;
400 libMesh::DenseVector<libMesh::Real> X(2); X(0) = X(1) = 0.0;
401 libMesh::DenseVector<libMesh::Real> dX(2);
404 libMesh::DenseVector<libMesh::Real> R(2);
405 libMesh::DenseMatrix<libMesh::Real> dRdX(2,2);
406 libMesh::DenseMatrix<libMesh::Real> dRdXinv(2,2);
411 libMesh::Real R0 = R.l2_norm();
414 while ( (R.l2_norm() > ATOL) && (R.l2_norm()/R0 > RTOL) && (iter<ITER_MAX) )
417 jac(X, R, res, dRdX);
420 det = dRdX(0,0)*dRdX(1,1) - dRdX(0,1)*dRdX(1,0);
423 if(std::abs(det)<=1e-20)
425 std::cout <<
"WARNING: about to divide by zero." << std::endl;
428 dRdXinv(0,0) = dRdX(1,1)/det;
429 dRdXinv(0,1) = -dRdX(0,1)/det;
430 dRdXinv(1,0) = -dRdX(1,0)/det;
431 dRdXinv(1,1) = dRdX(0,0)/det;
434 dX(0) = -( dRdXinv(0,0)*R(0) + dRdXinv(0,1)*R(1) );
435 dX(1) = -( dRdXinv(1,0)*R(0) + dRdXinv(1,1)*R(1) );
461 std::cout <<
"WARNING: Failed to converge distance function iteration. Assuming that min is outside this face. Use closest distance to edges to proceed." << std::endl;
469 if ( std::abs(X(0)) > 1.0 || std::abs(X(1)) > 1.0 )
474 libMesh::Real dtmp = std::numeric_limits<libMesh::Real>::infinity();
482 for (
unsigned int iedge=0; iedge<4; iedge++ ){
487 if (iedge==0) { line.p0 = belem->point(0); line.p1 = belem->point(1); }
488 else if (iedge==1) { line.p0 = belem->point(1); line.p1 = belem->point(2); }
489 else if (iedge==2) { line.p0 = belem->point(2); line.p1 = belem->point(3); }
490 else { line.p0 = belem->point(3); line.p1 = belem->point(0); }
492 DistanceToSegment<3>(*node, line, dtmp);
494 if( dtmp < dedge ) dedge = dtmp;
503 libMesh::UniquePtr<libMesh::FEBase> fe( libMesh::FEBase::build(2, libMesh::FEType(libMesh::FIRST, libMesh::LAGRANGE)) );
505 std::vector<libMesh::Point> xi(1);
510 const std::vector<std::vector<libMesh::Real> > &basis = fe->get_phi();
513 fe->reinit(belem, &xi);
516 libMesh::DenseVector<libMesh::Real> xx(3);
519 for (
unsigned int inode=0; inode<belem->n_nodes(); ++inode)
520 for (
unsigned int idim=0; idim<3; ++idim)
521 xx(idim) += belem->point(inode)(idim) * basis[inode][0];
525 for (
unsigned int ii=0; ii<3; ii++ ) {
526 dedge += (xnode[ii] - xx(ii))*(xnode[ii] - xx(ii));
535 std::cout <<
"My type is " << belem->type() << std::endl;
536 libmesh_not_implemented();
539 if ( dedge < distance ) distance = dedge;
563 libMesh::System& system =
_equation_systems.get_system<libMesh::System>(
"distance_function");
564 const unsigned int sys_num = system.number();
575 libMesh::MeshSerializer serialize(const_cast<libMesh::UnstructuredMesh&>(
_boundary_mesh));
578 libMesh::MeshBase::const_node_iterator node_it = mesh.local_nodes_begin();
579 const libMesh::MeshBase::const_node_iterator node_end = mesh.local_nodes_end();
588 for ( ; node_it != node_end; ++node_it) {
591 const libMesh::Node* node = *node_it;
597 const unsigned int dof = node->dof_number(sys_num,0,0);
598 system.solution->set (dof, distance);
604 system.solution->close();
607 std::cout <<
"Distance Function computed." << std::endl;
614 libMesh::UniquePtr< libMesh::DenseVector<libMesh::Real> >
617 libmesh_assert( elem != NULL );
618 libmesh_assert( qpts.size() > 0 );
621 const std::vector<std::vector<libMesh::Real> > &phi =
_dist_fe->get_phi();
627 const unsigned int n_dofs = phi.size();
630 const unsigned int n_pts = qpts.size();
633 libMesh::UniquePtr< DenseVector<libMesh::Real> > ap(
new libMesh::DenseVector<libMesh::Real>(qpts.size()) );
637 libMesh::System& sys =
_equation_systems.get_system<libMesh::System>(
"distance_function");
638 const libMesh::DofMap& dof_map = sys.get_dof_map();
640 std::vector<unsigned int> dof_ind;
641 dof_map.dof_indices(elem, dof_ind);
643 libMesh::DenseVector<libMesh::Real> nodal_dist;
644 nodal_dist.resize(n_dofs);
646 dof_map.extract_local_vector( *(sys.current_local_solution), dof_ind, nodal_dist);
648 for (
unsigned int idof=0; idof<n_dofs; idof++ ) {
649 for (
unsigned int iqpt=0; iqpt<n_pts; iqpt++ ) {
650 (*ap)(iqpt) += nodal_dist(idof) * phi[idof][iqpt];
libMesh::UniquePtr< libMesh::FEBase > _dist_fe
libMesh::EquationSystems & _equation_systems
DistanceFunction(libMesh::EquationSystems &es_in, const libMesh::UnstructuredMesh &bm_in)
libMesh::UniquePtr< libMesh::DenseVector< libMesh::Real > > interpolate(const libMesh::Elem *elem, const std::vector< libMesh::Point > &qts) const
virtual void initialize()
libMesh::Real node_to_boundary(const libMesh::Node *node)
const libMesh::UnstructuredMesh & _boundary_mesh