$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 // C++ includes 00021 00022 00023 // Local Includes 00024 #include "libmesh/mesh_function.h" 00025 #include "libmesh/dense_vector.h" 00026 #include "libmesh/equation_systems.h" 00027 #include "libmesh/numeric_vector.h" 00028 #include "libmesh/dof_map.h" 00029 #include "libmesh/point_locator_tree.h" 00030 #include "libmesh/fe_base.h" 00031 #include "libmesh/fe_interface.h" 00032 #include "libmesh/fe_compute_data.h" 00033 #include "libmesh/mesh_base.h" 00034 #include "libmesh/point.h" 00035 00036 namespace libMesh 00037 { 00038 00039 00040 //------------------------------------------------------------------ 00041 // MeshFunction methods 00042 MeshFunction::MeshFunction (const EquationSystems& eqn_systems, 00043 const NumericVector<Number>& vec, 00044 const DofMap& dof_map, 00045 const std::vector<unsigned int>& vars, 00046 const FunctionBase<Number>* master) : 00047 FunctionBase<Number> (master), 00048 ParallelObject (eqn_systems), 00049 _eqn_systems (eqn_systems), 00050 _vector (vec), 00051 _dof_map (dof_map), 00052 _system_vars (vars), 00053 _point_locator (NULL), 00054 _out_of_mesh_mode (false), 00055 _out_of_mesh_value () 00056 { 00057 } 00058 00059 00060 00061 MeshFunction::MeshFunction (const EquationSystems& eqn_systems, 00062 const NumericVector<Number>& vec, 00063 const DofMap& dof_map, 00064 const unsigned int var, 00065 const FunctionBase<Number>* master) : 00066 FunctionBase<Number> (master), 00067 ParallelObject (eqn_systems), 00068 _eqn_systems (eqn_systems), 00069 _vector (vec), 00070 _dof_map (dof_map), 00071 _system_vars (1,var), 00072 _point_locator (NULL), 00073 _out_of_mesh_mode (false), 00074 _out_of_mesh_value () 00075 { 00076 // std::vector<unsigned int> buf (1); 00077 // buf[0] = var; 00078 // _system_vars (buf); 00079 } 00080 00081 00082 00083 00084 00085 00086 00087 MeshFunction::~MeshFunction () 00088 { 00089 // only delete the point locator when we are the master 00090 if (this->_master == NULL) 00091 delete this->_point_locator; 00092 } 00093 00094 00095 00096 00097 void MeshFunction::init (const Trees::BuildType /*point_locator_build_type*/) 00098 { 00099 // are indices of the desired variable(s) provided? 00100 libmesh_assert_greater (this->_system_vars.size(), 0); 00101 00102 // Don't do twice... 00103 if (this->_initialized) 00104 { 00105 libmesh_assert(this->_point_locator); 00106 return; 00107 } 00108 00109 /* 00110 * set up the PointLocator: either someone else 00111 * is the master (go and get the address of his 00112 * point locator) or this object is the master 00113 * (build the point locator on our own). 00114 */ 00115 if (this->_master != NULL) 00116 { 00117 // we aren't the master 00118 const MeshFunction* master = 00119 cast_ptr<const MeshFunction*>(this->_master); 00120 00121 if (master->_point_locator == NULL) 00122 libmesh_error_msg("ERROR: When the master-servant concept is used, the master has to be initialized first!"); 00123 00124 else 00125 { 00126 this->_point_locator = master->_point_locator; 00127 } 00128 } 00129 else 00130 { 00131 // we are the master: build the point locator 00132 00133 // constant reference to the other mesh 00134 const MeshBase& mesh = this->_eqn_systems.get_mesh(); 00135 00136 // build the point locator. Only \p TREE version available 00137 //UniquePtr<PointLocatorBase> ap (PointLocatorBase::build (TREE, mesh)); 00138 //this->_point_locator = ap.release(); 00139 // this->_point_locator = new PointLocatorTree (mesh, point_locator_build_type); 00140 this->_point_locator = mesh.sub_point_locator().release(); 00141 00142 // Point locator no longer needs to be initialized. 00143 // this->_point_locator->init(); 00144 } 00145 00146 00147 // ready for use 00148 this->_initialized = true; 00149 } 00150 00151 00152 void 00153 MeshFunction::clear () 00154 { 00155 // only delete the point locator when we are the master 00156 if ((this->_point_locator != NULL) && (this->_master == NULL)) 00157 { 00158 delete this->_point_locator; 00159 this->_point_locator = NULL; 00160 } 00161 this->_initialized = false; 00162 } 00163 00164 00165 00166 UniquePtr<FunctionBase<Number> > MeshFunction::clone () const 00167 { 00168 return UniquePtr<FunctionBase<Number> > 00169 (new MeshFunction 00170 (_eqn_systems, _vector, _dof_map, _system_vars, this)); 00171 } 00172 00173 00174 00175 Number MeshFunction::operator() (const Point& p, 00176 const Real time) 00177 { 00178 libmesh_assert (this->initialized()); 00179 00180 DenseVector<Number> buf (1); 00181 this->operator() (p, time, buf); 00182 return buf(0); 00183 } 00184 00185 00186 00187 Gradient MeshFunction::gradient (const Point& p, 00188 const Real time) 00189 { 00190 libmesh_assert (this->initialized()); 00191 00192 std::vector<Gradient> buf (1); 00193 this->gradient(p, time, buf); 00194 return buf[0]; 00195 } 00196 00197 00198 00199 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00200 Tensor MeshFunction::hessian (const Point& p, 00201 const Real time) 00202 { 00203 libmesh_assert (this->initialized()); 00204 00205 std::vector<Tensor> buf (1); 00206 this->hessian(p, time, buf); 00207 return buf[0]; 00208 } 00209 #endif 00210 00211 void MeshFunction::operator() (const Point& p, 00212 const Real time, 00213 DenseVector<Number>& output) 00214 { 00215 this->operator() (p,time,output,NULL); 00216 } 00217 00218 void MeshFunction::operator() (const Point& p, 00219 const Real, 00220 DenseVector<Number>& output, 00221 const std::set<subdomain_id_type>* subdomain_ids) 00222 { 00223 libmesh_assert (this->initialized()); 00224 00225 const Elem* element = this->find_element(p,subdomain_ids); 00226 00227 if (!element) 00228 { 00229 output = _out_of_mesh_value; 00230 } 00231 else 00232 { 00233 // resize the output vector to the number of output values 00234 // that the user told us 00235 output.resize (cast_int<unsigned int> 00236 (this->_system_vars.size())); 00237 00238 00239 { 00240 const unsigned int dim = element->dim(); 00241 00242 00243 /* 00244 * Get local coordinates to feed these into compute_data(). 00245 * Note that the fe_type can safely be used from the 0-variable, 00246 * since the inverse mapping is the same for all FEFamilies 00247 */ 00248 const Point mapped_point (FEInterface::inverse_map (dim, 00249 this->_dof_map.variable_type(0), 00250 element, 00251 p)); 00252 00253 00254 // loop over all vars 00255 for (unsigned int index=0; index < this->_system_vars.size(); index++) 00256 { 00257 /* 00258 * the data for this variable 00259 */ 00260 const unsigned int var = _system_vars[index]; 00261 const FEType& fe_type = this->_dof_map.variable_type(var); 00262 00267 { 00268 FEComputeData data (this->_eqn_systems, mapped_point); 00269 00270 FEInterface::compute_data (dim, fe_type, element, data); 00271 00272 // where the solution values for the var-th variable are stored 00273 std::vector<dof_id_type> dof_indices; 00274 this->_dof_map.dof_indices (element, dof_indices, var); 00275 00276 // interpolate the solution 00277 { 00278 Number value = 0.; 00279 00280 for (unsigned int i=0; i<dof_indices.size(); i++) 00281 value += this->_vector(dof_indices[i]) * data.shape[i]; 00282 00283 output(index) = value; 00284 } 00285 00286 } 00287 00288 // next variable 00289 } 00290 } 00291 } 00292 00293 // all done 00294 return; 00295 } 00296 00297 00298 00299 void MeshFunction::gradient (const Point& p, 00300 const Real, 00301 std::vector<Gradient>& output, 00302 const std::set<subdomain_id_type>* subdomain_ids) 00303 { 00304 libmesh_assert (this->initialized()); 00305 00306 const Elem* element = this->find_element(p,subdomain_ids); 00307 00308 if (!element) 00309 { 00310 output.resize(0); 00311 } 00312 else 00313 { 00314 // resize the output vector to the number of output values 00315 // that the user told us 00316 output.resize (this->_system_vars.size()); 00317 00318 00319 { 00320 const unsigned int dim = element->dim(); 00321 00322 00323 /* 00324 * Get local coordinates to feed these into compute_data(). 00325 * Note that the fe_type can safely be used from the 0-variable, 00326 * since the inverse mapping is the same for all FEFamilies 00327 */ 00328 const Point mapped_point (FEInterface::inverse_map (dim, 00329 this->_dof_map.variable_type(0), 00330 element, 00331 p)); 00332 00333 std::vector<Point> point_list (1, mapped_point); 00334 00335 // loop over all vars 00336 for (unsigned int index=0; index < this->_system_vars.size(); index++) 00337 { 00338 /* 00339 * the data for this variable 00340 */ 00341 const unsigned int var = _system_vars[index]; 00342 const FEType& fe_type = this->_dof_map.variable_type(var); 00343 00344 UniquePtr<FEBase> point_fe (FEBase::build(dim, fe_type)); 00345 const std::vector<std::vector<RealGradient> >& dphi = point_fe->get_dphi(); 00346 point_fe->reinit(element, &point_list); 00347 00348 // where the solution values for the var-th variable are stored 00349 std::vector<dof_id_type> dof_indices; 00350 this->_dof_map.dof_indices (element, dof_indices, var); 00351 00352 // interpolate the solution 00353 Gradient grad(0.); 00354 00355 for (unsigned int i=0; i<dof_indices.size(); i++) 00356 grad.add_scaled(dphi[i][0], this->_vector(dof_indices[i])); 00357 00358 output[index] = grad; 00359 } 00360 } 00361 } 00362 00363 // all done 00364 return; 00365 } 00366 00367 00368 00369 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00370 void MeshFunction::hessian (const Point& p, 00371 const Real, 00372 std::vector<Tensor>& output, 00373 const std::set<subdomain_id_type>* subdomain_ids) 00374 { 00375 libmesh_assert (this->initialized()); 00376 00377 const Elem* element = this->find_element(p,subdomain_ids); 00378 00379 if (!element) 00380 { 00381 output.resize(0); 00382 } 00383 else 00384 { 00385 // resize the output vector to the number of output values 00386 // that the user told us 00387 output.resize (this->_system_vars.size()); 00388 00389 00390 { 00391 const unsigned int dim = element->dim(); 00392 00393 00394 /* 00395 * Get local coordinates to feed these into compute_data(). 00396 * Note that the fe_type can safely be used from the 0-variable, 00397 * since the inverse mapping is the same for all FEFamilies 00398 */ 00399 const Point mapped_point (FEInterface::inverse_map (dim, 00400 this->_dof_map.variable_type(0), 00401 element, 00402 p)); 00403 00404 std::vector<Point> point_list (1, mapped_point); 00405 00406 // loop over all vars 00407 for (unsigned int index=0; index < this->_system_vars.size(); index++) 00408 { 00409 /* 00410 * the data for this variable 00411 */ 00412 const unsigned int var = _system_vars[index]; 00413 const FEType& fe_type = this->_dof_map.variable_type(var); 00414 00415 UniquePtr<FEBase> point_fe (FEBase::build(dim, fe_type)); 00416 const std::vector<std::vector<RealTensor> >& d2phi = 00417 point_fe->get_d2phi(); 00418 point_fe->reinit(element, &point_list); 00419 00420 // where the solution values for the var-th variable are stored 00421 std::vector<dof_id_type> dof_indices; 00422 this->_dof_map.dof_indices (element, dof_indices, var); 00423 00424 // interpolate the solution 00425 Tensor hess; 00426 00427 for (unsigned int i=0; i<dof_indices.size(); i++) 00428 hess.add_scaled(d2phi[i][0], this->_vector(dof_indices[i])); 00429 00430 output[index] = hess; 00431 } 00432 } 00433 } 00434 00435 // all done 00436 return; 00437 } 00438 #endif 00439 00440 const Elem* MeshFunction::find_element( const Point& p, 00441 const std::set<subdomain_id_type>* subdomain_ids ) const 00442 { 00443 /* Ensure that in the case of a master mesh function, the 00444 out-of-mesh mode is enabled either for both or for none. This is 00445 important because the out-of-mesh mode is also communicated to 00446 the point locator. Since this is time consuming, enable it only 00447 in debug mode. */ 00448 #ifdef DEBUG 00449 if (this->_master != NULL) 00450 { 00451 const MeshFunction* master = 00452 cast_ptr<const MeshFunction*>(this->_master); 00453 if(_out_of_mesh_mode!=master->_out_of_mesh_mode) 00454 libmesh_error_msg("ERROR: If you use out-of-mesh-mode in connection with master mesh " \ 00455 << "functions, you must enable out-of-mesh mode for both the master and the slave mesh function."); 00456 } 00457 #endif 00458 00459 // locate the point in the other mesh 00460 const Elem* element = this->_point_locator->operator()(p,subdomain_ids); 00461 00462 // If we have an element, but it's not a local element, then we 00463 // either need to have a serialized vector or we need to find a 00464 // local element sharing the same point. 00465 if (element && 00466 (element->processor_id() != this->processor_id()) && 00467 _vector.type() != SERIAL) 00468 { 00469 // look for a local element containing the point 00470 std::set<const Elem*> point_neighbors; 00471 element->find_point_neighbors(p, point_neighbors); 00472 element = NULL; 00473 std::set<const Elem*>::const_iterator it = point_neighbors.begin(); 00474 const std::set<const Elem*>::const_iterator end = point_neighbors.end(); 00475 for (; it != end; ++it) 00476 { 00477 const Elem* elem = *it; 00478 if (elem->processor_id() == this->processor_id()) 00479 { 00480 element = elem; 00481 break; 00482 } 00483 } 00484 } 00485 00486 return element; 00487 } 00488 00489 const PointLocatorBase& MeshFunction::get_point_locator (void) const 00490 { 00491 libmesh_assert (this->initialized()); 00492 return *_point_locator; 00493 } 00494 00495 void MeshFunction::enable_out_of_mesh_mode(const DenseVector<Number>& value) 00496 { 00497 libmesh_assert (this->initialized()); 00498 _point_locator->enable_out_of_mesh_mode(); 00499 _out_of_mesh_mode = true; 00500 _out_of_mesh_value = value; 00501 } 00502 00503 void MeshFunction::enable_out_of_mesh_mode(const Number& value) 00504 { 00505 DenseVector<Number> v(1); 00506 v(0) = value; 00507 this->enable_out_of_mesh_mode(v); 00508 } 00509 00510 void MeshFunction::disable_out_of_mesh_mode(void) 00511 { 00512 libmesh_assert (this->initialized()); 00513 _point_locator->disable_out_of_mesh_mode(); 00514 _out_of_mesh_mode = false; 00515 } 00516 00517 void MeshFunction::set_point_locator_tolerance(Real tol) 00518 { 00519 // We need to enable out_of_mesh mode in the point_locator 00520 // in order for the point locator tolerance to be used. 00521 _point_locator->enable_out_of_mesh_mode(); 00522 00523 // Set the tolerance 00524 _point_locator->set_close_to_point_tol(tol); 00525 } 00526 00527 void MeshFunction::unset_point_locator_tolerance() 00528 { 00529 _point_locator->unset_close_to_point_tol(); 00530 } 00531 00532 } // namespace libMesh