system.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 #include <sstream> // for std::ostringstream 00022 00023 00024 // Local includes 00025 #include "libmesh/dof_map.h" 00026 #include "libmesh/equation_systems.h" 00027 #include "libmesh/libmesh_logging.h" 00028 #include "libmesh/mesh_base.h" 00029 #include "libmesh/numeric_vector.h" 00030 #include "libmesh/parameter_vector.h" 00031 #include "libmesh/point.h" // For point_value 00032 #include "libmesh/point_locator_base.h" // For point_value 00033 #include "libmesh/qoi_set.h" 00034 #include "libmesh/string_to_enum.h" 00035 #include "libmesh/system.h" 00036 #include "libmesh/system_norm.h" 00037 #include "libmesh/utility.h" 00038 00039 // includes for calculate_norm, point_* 00040 #include "libmesh/fe_base.h" 00041 #include "libmesh/fe_interface.h" 00042 #include "libmesh/parallel.h" 00043 #include "libmesh/parallel_algebra.h" 00044 #include "libmesh/quadrature.h" 00045 #include "libmesh/tensor_value.h" 00046 #include "libmesh/vector_value.h" 00047 #include "libmesh/tensor_tools.h" 00048 00049 namespace libMesh 00050 { 00051 00052 00053 // ------------------------------------------------------------ 00054 // System implementation 00055 System::System (EquationSystems& es, 00056 const std::string& name_in, 00057 const unsigned int number_in) : 00058 00059 assemble_before_solve (true), 00060 use_fixed_solution (false), 00061 extra_quadrature_order (0), 00062 solution (NumericVector<Number>::build()), 00063 current_local_solution (NumericVector<Number>::build()), 00064 time (0.), 00065 qoi (0), 00066 _init_system_function (NULL), 00067 _init_system_object (NULL), 00068 _assemble_system_function (NULL), 00069 _assemble_system_object (NULL), 00070 _constrain_system_function (NULL), 00071 _constrain_system_object (NULL), 00072 _qoi_evaluate_function (NULL), 00073 _qoi_evaluate_object (NULL), 00074 _qoi_evaluate_derivative_function (NULL), 00075 _qoi_evaluate_derivative_object (NULL), 00076 _dof_map (new DofMap(number_in)), 00077 _equation_systems (es), 00078 _mesh (es.get_mesh()), 00079 _sys_name (name_in), 00080 _sys_number (number_in), 00081 _active (true), 00082 _solution_projection (true), 00083 _basic_system_only (false), 00084 _can_add_vectors (true), 00085 _identify_variable_groups (true), 00086 _additional_data_written (false) 00087 { 00088 } 00089 00090 00091 00092 // No copy construction of System objects! 00093 System::System (const System& other) : 00094 ReferenceCountedObject<System>(), 00095 _equation_systems(other._equation_systems), 00096 _mesh(other._mesh), 00097 _sys_number(other._sys_number) 00098 { 00099 libmesh_error(); 00100 } 00101 00102 00103 00104 System& System::operator= (const System&) 00105 { 00106 libmesh_error(); 00107 } 00108 00109 00110 System::~System () 00111 { 00112 // Null-out the function pointers. Since this 00113 // class is getting destructed it is pointless, 00114 // but a good habit. 00115 _init_system_function = 00116 _assemble_system_function = 00117 _constrain_system_function = NULL; 00118 00119 _qoi_evaluate_function = 00120 _qoi_evaluate_derivative_function = NULL; 00121 00122 // NULL-out user-provided objects. 00123 _init_system_object = NULL; 00124 _assemble_system_object = NULL; 00125 _constrain_system_object = NULL; 00126 _qoi_evaluate_object = NULL; 00127 _qoi_evaluate_derivative_object = NULL; 00128 00129 // Clear data 00130 this->clear (); 00131 00132 libmesh_assert (!libMesh::closed()); 00133 } 00134 00135 00136 00137 dof_id_type System::n_dofs() const 00138 { 00139 return _dof_map->n_dofs(); 00140 } 00141 00142 00143 00144 dof_id_type System::n_constrained_dofs() const 00145 { 00146 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00147 00148 return _dof_map->n_constrained_dofs(); 00149 00150 #else 00151 00152 return 0; 00153 00154 #endif 00155 } 00156 00157 00158 00159 dof_id_type System::n_local_constrained_dofs() const 00160 { 00161 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00162 00163 return _dof_map->n_local_constrained_dofs(); 00164 00165 #else 00166 00167 return 0; 00168 00169 #endif 00170 } 00171 00172 00173 00174 dof_id_type System::n_local_dofs() const 00175 { 00176 return _dof_map->n_dofs_on_processor (libMesh::processor_id()); 00177 } 00178 00179 00180 00181 Number System::current_solution (const dof_id_type global_dof_number) const 00182 { 00183 // Check the sizes 00184 libmesh_assert_less (global_dof_number, _dof_map->n_dofs()); 00185 libmesh_assert_less (global_dof_number, current_local_solution->size()); 00186 00187 return (*current_local_solution)(global_dof_number); 00188 } 00189 00190 00191 00192 void System::clear () 00193 { 00194 _variables.clear(); 00195 00196 _variable_numbers.clear(); 00197 00198 _dof_map->clear (); 00199 00200 solution->clear (); 00201 00202 current_local_solution->clear (); 00203 00204 // clear any user-added vectors 00205 { 00206 for (vectors_iterator pos = _vectors.begin(); pos != _vectors.end(); ++pos) 00207 { 00208 pos->second->clear (); 00209 delete pos->second; 00210 pos->second = NULL; 00211 } 00212 00213 _vectors.clear(); 00214 _vector_projections.clear(); 00215 _vector_types.clear(); 00216 _can_add_vectors = true; 00217 } 00218 00219 } 00220 00221 00222 00223 void System::init () 00224 { 00225 // First initialize any required data: 00226 // either only the basic System data 00227 if (_basic_system_only) 00228 System::init_data(); 00229 // or all the derived class' data too 00230 else 00231 this->init_data(); 00232 00233 // If no variables have been added to this system 00234 // don't do anything 00235 if(!this->n_vars()) 00236 return; 00237 00238 // Then call the user-provided intialization function 00239 this->user_initialization(); 00240 } 00241 00242 00243 00244 void System::init_data () 00245 { 00246 MeshBase &mesh = this->get_mesh(); 00247 00248 // Add all variable groups to our underlying DofMap 00249 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 00250 _dof_map->add_variable_group(this->variable_group(vg)); 00251 00252 // Distribute the degrees of freedom on the mesh 00253 _dof_map->distribute_dofs (mesh); 00254 00255 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00256 00257 // Recreate any hanging node constraints 00258 _dof_map->create_dof_constraints(mesh); 00259 00260 // Apply any user-defined constraints 00261 this->user_constrain(); 00262 00263 // Expand any recursive constraints 00264 _dof_map->process_constraints(mesh); 00265 00266 #endif 00267 00268 // And clean up the send_list before we first use it 00269 _dof_map->prepare_send_list(); 00270 00271 // Resize the solution conformal to the current mesh 00272 solution->init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL); 00273 00274 // Resize the current_local_solution for the current mesh 00275 #ifdef LIBMESH_ENABLE_GHOSTED 00276 current_local_solution->init (this->n_dofs(), this->n_local_dofs(), 00277 _dof_map->get_send_list(), false, 00278 GHOSTED); 00279 #else 00280 current_local_solution->init (this->n_dofs(), false, SERIAL); 00281 #endif 00282 00283 // from now on, no chance to add additional vectors 00284 _can_add_vectors = false; 00285 00286 // initialize & zero other vectors, if necessary 00287 for (vectors_iterator pos = _vectors.begin(); pos != _vectors.end(); ++pos) 00288 { 00289 ParallelType type = _vector_types[pos->first]; 00290 00291 if (type == GHOSTED) 00292 { 00293 #ifdef LIBMESH_ENABLE_GHOSTED 00294 pos->second->init (this->n_dofs(), this->n_local_dofs(), 00295 _dof_map->get_send_list(), false, 00296 GHOSTED); 00297 #else 00298 libMesh::err << "Cannot initialize ghosted vectors when they are not enabled." << std::endl; 00299 libmesh_error(); 00300 #endif 00301 } 00302 else if (type == SERIAL) 00303 { 00304 pos->second->init (this->n_dofs(), false, type); 00305 } 00306 else 00307 { 00308 libmesh_assert_equal_to(type, PARALLEL); 00309 pos->second->init (this->n_dofs(), this->n_local_dofs(), false, type); 00310 } 00311 } 00312 } 00313 00314 00315 00316 void System::restrict_vectors () 00317 { 00318 #ifdef LIBMESH_ENABLE_AMR 00319 // Restrict the _vectors on the coarsened cells 00320 for (vectors_iterator pos = _vectors.begin(); pos != _vectors.end(); ++pos) 00321 { 00322 NumericVector<Number>* v = pos->second; 00323 00324 if (_vector_projections[pos->first]) 00325 this->project_vector (*v); 00326 else 00327 { 00328 ParallelType type = _vector_types[pos->first]; 00329 00330 if(type == GHOSTED) 00331 { 00332 #ifdef LIBMESH_ENABLE_GHOSTED 00333 pos->second->init (this->n_dofs(), this->n_local_dofs(), 00334 _dof_map->get_send_list(), false, 00335 GHOSTED); 00336 #else 00337 libMesh::err << "Cannot initialize ghosted vectors when they are not enabled." << std::endl; 00338 libmesh_error(); 00339 #endif 00340 } 00341 else 00342 pos->second->init (this->n_dofs(), this->n_local_dofs(), false, type); 00343 } 00344 } 00345 00346 const std::vector<dof_id_type>& send_list = _dof_map->get_send_list (); 00347 00348 // Restrict the solution on the coarsened cells 00349 if (_solution_projection) 00350 this->project_vector (*solution); 00351 00352 #ifdef LIBMESH_ENABLE_GHOSTED 00353 current_local_solution->init(this->n_dofs(), 00354 this->n_local_dofs(), send_list, 00355 false, GHOSTED); 00356 #else 00357 current_local_solution->init(this->n_dofs()); 00358 #endif 00359 00360 if (_solution_projection) 00361 solution->localize (*current_local_solution, send_list); 00362 00363 #endif // LIBMESH_ENABLE_AMR 00364 } 00365 00366 00367 00368 void System::prolong_vectors () 00369 { 00370 #ifdef LIBMESH_ENABLE_AMR 00371 // Currently project_vector handles both restriction and prolongation 00372 this->restrict_vectors(); 00373 #endif 00374 } 00375 00376 00377 00378 void System::reinit () 00379 { 00380 #ifdef LIBMESH_ENABLE_AMR 00381 //If no variables have been added to this system 00382 //don't do anything 00383 if(!this->n_vars()) 00384 return; 00385 00386 // Constraints get handled in EquationSystems::reinit now 00387 // _dof_map->create_dof_constraints(this->get_mesh()); 00388 00389 // Update the solution based on the projected 00390 // current_local_solution. 00391 solution->init (this->n_dofs(), this->n_local_dofs(), true, PARALLEL); 00392 00393 libmesh_assert_equal_to (solution->size(), current_local_solution->size()); 00394 // Not true with ghosted vectors 00395 // libmesh_assert_equal_to (solution->size(), current_local_solution->local_size()); 00396 00397 const dof_id_type first_local_dof = solution->first_local_index(); 00398 const dof_id_type local_size = solution->local_size(); 00399 00400 for (dof_id_type i=0; i<local_size; i++) 00401 solution->set(i+first_local_dof, 00402 (*current_local_solution)(i+first_local_dof)); 00403 00404 solution->close(); 00405 #endif 00406 } 00407 00408 00409 00410 void System::update () 00411 { 00412 libmesh_assert(solution->closed()); 00413 00414 const std::vector<dof_id_type>& send_list = _dof_map->get_send_list (); 00415 00416 // Check sizes 00417 libmesh_assert_equal_to (current_local_solution->size(), solution->size()); 00418 // More processors than elements => empty send_list 00419 // libmesh_assert (!send_list.empty()); 00420 libmesh_assert_less_equal (send_list.size(), solution->size()); 00421 00422 // Create current_local_solution from solution. This will 00423 // put a local copy of solution into current_local_solution. 00424 // Only the necessary values (specified by the send_list) 00425 // are copied to minimize communication 00426 solution->localize (*current_local_solution, send_list); 00427 } 00428 00429 00430 00431 void System::re_update () 00432 { 00433 parallel_only(); 00434 00435 // If this system is empty... don't do anything! 00436 if(!this->n_vars()) 00437 return; 00438 00439 const std::vector<dof_id_type>& send_list = this->get_dof_map().get_send_list (); 00440 00441 // Check sizes 00442 libmesh_assert_equal_to (current_local_solution->size(), solution->size()); 00443 // Not true with ghosted vectors 00444 // libmesh_assert_equal_to (current_local_solution->local_size(), solution->size()); 00445 // libmesh_assert (!send_list.empty()); 00446 libmesh_assert_less_equal (send_list.size(), solution->size()); 00447 00448 // Create current_local_solution from solution. This will 00449 // put a local copy of solution into current_local_solution. 00450 solution->localize (*current_local_solution, send_list); 00451 } 00452 00453 00454 00455 void System::restrict_solve_to (const SystemSubset* subset, 00456 const SubsetSolveMode /*subset_solve_mode*/) 00457 { 00458 if(subset!=NULL) 00459 { 00460 libmesh_not_implemented(); 00461 } 00462 } 00463 00464 00465 00466 void System::assemble () 00467 { 00468 // Log how long the user's assembly code takes 00469 START_LOG("assemble()", "System"); 00470 00471 // Call the user-specified assembly function 00472 this->user_assembly(); 00473 00474 // Stop logging the user code 00475 STOP_LOG("assemble()", "System"); 00476 } 00477 00478 00479 00480 void System::assemble_qoi (const QoISet& qoi_indices) 00481 { 00482 // Log how long the user's assembly code takes 00483 START_LOG("assemble_qoi()", "System"); 00484 00485 // Call the user-specified quantity of interest function 00486 this->user_QOI(qoi_indices); 00487 00488 // Stop logging the user code 00489 STOP_LOG("assemble_qoi()", "System"); 00490 } 00491 00492 00493 00494 void System::assemble_qoi_derivative (const QoISet& qoi_indices) 00495 { 00496 // Log how long the user's assembly code takes 00497 START_LOG("assemble_qoi_derivative()", "System"); 00498 00499 // Call the user-specified quantity of interest function 00500 this->user_QOI_derivative(qoi_indices); 00501 00502 // Stop logging the user code 00503 STOP_LOG("assemble_qoi_derivative()", "System"); 00504 } 00505 00506 00507 00508 void System::qoi_parameter_sensitivity 00509 (const QoISet& qoi_indices, 00510 const ParameterVector& parameters, 00511 SensitivityData& sensitivities) 00512 { 00513 // Forward sensitivities are more efficient for Nq > Np 00514 if (qoi_indices.size(*this) > parameters.size()) 00515 forward_qoi_parameter_sensitivity(qoi_indices, parameters, sensitivities); 00516 // Adjoint sensitivities are more efficient for Np > Nq, 00517 // and an adjoint may be more reusable than a forward 00518 // solution sensitivity in the Np == Nq case. 00519 else 00520 adjoint_qoi_parameter_sensitivity(qoi_indices, parameters, sensitivities); 00521 } 00522 00523 00524 00525 bool System::compare (const System& other_system, 00526 const Real threshold, 00527 const bool verbose) const 00528 { 00529 // we do not care for matrices, but for vectors 00530 libmesh_assert (!_can_add_vectors); 00531 libmesh_assert (!other_system._can_add_vectors); 00532 00533 if (verbose) 00534 { 00535 libMesh::out << " Systems \"" << _sys_name << "\"" << std::endl; 00536 libMesh::out << " comparing matrices not supported." << std::endl; 00537 libMesh::out << " comparing names..."; 00538 } 00539 00540 // compare the name: 0 means identical 00541 const int name_result = _sys_name.compare(other_system.name()); 00542 if (verbose) 00543 { 00544 if (name_result == 0) 00545 libMesh::out << " identical." << std::endl; 00546 else 00547 libMesh::out << " names not identical." << std::endl; 00548 libMesh::out << " comparing solution vector..."; 00549 } 00550 00551 00552 // compare the solution: -1 means identical 00553 const int solu_result = solution->compare (*other_system.solution.get(), 00554 threshold); 00555 00556 if (verbose) 00557 { 00558 if (solu_result == -1) 00559 libMesh::out << " identical up to threshold." << std::endl; 00560 else 00561 libMesh::out << " first difference occured at index = " 00562 << solu_result << "." << std::endl; 00563 } 00564 00565 00566 // safety check, whether we handle at least the same number 00567 // of vectors 00568 std::vector<int> ov_result; 00569 00570 if (this->n_vectors() != other_system.n_vectors()) 00571 { 00572 if (verbose) 00573 { 00574 libMesh::out << " Fatal difference. This system handles " 00575 << this->n_vectors() << " add'l vectors," << std::endl 00576 << " while the other system handles " 00577 << other_system.n_vectors() 00578 << " add'l vectors." << std::endl 00579 << " Aborting comparison." << std::endl; 00580 } 00581 return false; 00582 } 00583 else if (this->n_vectors() == 0) 00584 { 00585 // there are no additional vectors... 00586 ov_result.clear (); 00587 } 00588 else 00589 { 00590 // compare other vectors 00591 for (const_vectors_iterator pos = _vectors.begin(); 00592 pos != _vectors.end(); ++pos) 00593 { 00594 if (verbose) 00595 libMesh::out << " comparing vector \"" 00596 << pos->first << "\" ..."; 00597 00598 // assume they have the same name 00599 const NumericVector<Number>& other_system_vector = 00600 other_system.get_vector(pos->first); 00601 00602 ov_result.push_back(pos->second->compare (other_system_vector, 00603 threshold)); 00604 00605 if (verbose) 00606 { 00607 if (ov_result[ov_result.size()-1] == -1) 00608 libMesh::out << " identical up to threshold." << std::endl; 00609 else 00610 libMesh::out << " first difference occured at" << std::endl 00611 << " index = " << ov_result[ov_result.size()-1] << "." << std::endl; 00612 } 00613 00614 } 00615 00616 } // finished comparing additional vectors 00617 00618 00619 bool overall_result; 00620 00621 // sum up the results 00622 if ((name_result==0) && (solu_result==-1)) 00623 { 00624 if (ov_result.size()==0) 00625 overall_result = true; 00626 else 00627 { 00628 bool ov_identical; 00629 unsigned int n = 0; 00630 do 00631 { 00632 ov_identical = (ov_result[n]==-1); 00633 n++; 00634 } 00635 while (ov_identical && n<ov_result.size()); 00636 overall_result = ov_identical; 00637 } 00638 } 00639 else 00640 overall_result = false; 00641 00642 if (verbose) 00643 { 00644 libMesh::out << " finished comparisons, "; 00645 if (overall_result) 00646 libMesh::out << "found no differences." << std::endl << std::endl; 00647 else 00648 libMesh::out << "found differences." << std::endl << std::endl; 00649 } 00650 00651 return overall_result; 00652 } 00653 00654 00655 00656 void System::update_global_solution (std::vector<Number>& global_soln) const 00657 { 00658 global_soln.resize (solution->size()); 00659 00660 solution->localize (global_soln); 00661 } 00662 00663 00664 00665 void System::update_global_solution (std::vector<Number>& global_soln, 00666 const unsigned int dest_proc) const 00667 { 00668 global_soln.resize (solution->size()); 00669 00670 solution->localize_to_one (global_soln, dest_proc); 00671 } 00672 00673 00674 00675 NumericVector<Number> & System::add_vector (const std::string& vec_name, 00676 const bool projections, 00677 const ParallelType type) 00678 { 00679 // Return the vector if it is already there. 00680 if (this->have_vector(vec_name)) 00681 return *(_vectors[vec_name]); 00682 00683 // Otherwise build the vector 00684 NumericVector<Number>* buf = NumericVector<Number>::build().release(); 00685 _vectors.insert (std::make_pair (vec_name, buf)); 00686 _vector_projections.insert (std::make_pair (vec_name, projections)); 00687 00688 _vector_types.insert (std::make_pair (vec_name, type)); 00689 00690 // Initialize it if necessary 00691 if (!_can_add_vectors) 00692 { 00693 if(type == GHOSTED) 00694 { 00695 #ifdef LIBMESH_ENABLE_GHOSTED 00696 buf->init (this->n_dofs(), this->n_local_dofs(), 00697 _dof_map->get_send_list(), false, 00698 GHOSTED); 00699 #else 00700 std::cerr << "Cannot initialize ghosted vectors when they are not enabled." << std::endl; 00701 libmesh_error(); 00702 #endif 00703 } 00704 else 00705 buf->init (this->n_dofs(), this->n_local_dofs(), false, type); 00706 } 00707 00708 return *buf; 00709 } 00710 00711 void System::remove_vector (const std::string& vec_name) 00712 { 00713 //Return if the vector does not exist 00714 if ( !(this->have_vector(vec_name)) ) 00715 return; 00716 00717 _vectors[vec_name]->clear(); 00718 delete _vectors[vec_name]; 00719 _vectors[vec_name] = NULL; 00720 00721 _vectors.erase(vec_name); 00722 _vector_projections.erase(vec_name); 00723 _vector_types.erase(vec_name); 00724 } 00725 00726 const NumericVector<Number> * System::request_vector (const std::string& vec_name) const 00727 { 00728 const_vectors_iterator pos = _vectors.find(vec_name); 00729 00730 if (pos == _vectors.end()) 00731 return NULL; 00732 00733 return pos->second; 00734 } 00735 00736 00737 00738 NumericVector<Number> * System::request_vector (const std::string& vec_name) 00739 { 00740 vectors_iterator pos = _vectors.find(vec_name); 00741 00742 if (pos == _vectors.end()) 00743 return NULL; 00744 00745 return pos->second; 00746 } 00747 00748 00749 00750 const NumericVector<Number> * System::request_vector (const unsigned int vec_num) const 00751 { 00752 const_vectors_iterator v = vectors_begin(); 00753 const_vectors_iterator v_end = vectors_end(); 00754 unsigned int num = 0; 00755 while((num<vec_num) && (v!=v_end)) 00756 { 00757 num++; 00758 ++v; 00759 } 00760 if (v==v_end) 00761 return NULL; 00762 return v->second; 00763 } 00764 00765 00766 00767 NumericVector<Number> * System::request_vector (const unsigned int vec_num) 00768 { 00769 vectors_iterator v = vectors_begin(); 00770 vectors_iterator v_end = vectors_end(); 00771 unsigned int num = 0; 00772 while((num<vec_num) && (v!=v_end)) 00773 { 00774 num++; 00775 ++v; 00776 } 00777 if (v==v_end) 00778 return NULL; 00779 return v->second; 00780 } 00781 00782 00783 00784 const NumericVector<Number> & System::get_vector (const std::string& vec_name) const 00785 { 00786 // Make sure the vector exists 00787 const_vectors_iterator pos = _vectors.find(vec_name); 00788 00789 if (pos == _vectors.end()) 00790 { 00791 libMesh::err << "ERROR: vector " 00792 << vec_name 00793 << " does not exist in this system!" 00794 << std::endl; 00795 libmesh_error(); 00796 } 00797 00798 return *(pos->second); 00799 } 00800 00801 00802 00803 NumericVector<Number> & System::get_vector (const std::string& vec_name) 00804 { 00805 // Make sure the vector exists 00806 vectors_iterator pos = _vectors.find(vec_name); 00807 00808 if (pos == _vectors.end()) 00809 { 00810 libMesh::err << "ERROR: vector " 00811 << vec_name 00812 << " does not exist in this system!" 00813 << std::endl; 00814 libmesh_error(); 00815 } 00816 00817 return *(pos->second); 00818 } 00819 00820 00821 00822 const NumericVector<Number> & System::get_vector (const unsigned int vec_num) const 00823 { 00824 const_vectors_iterator v = vectors_begin(); 00825 const_vectors_iterator v_end = vectors_end(); 00826 unsigned int num = 0; 00827 while((num<vec_num) && (v!=v_end)) 00828 { 00829 num++; 00830 ++v; 00831 } 00832 libmesh_assert (v != v_end); 00833 return *(v->second); 00834 } 00835 00836 00837 00838 NumericVector<Number> & System::get_vector (const unsigned int vec_num) 00839 { 00840 vectors_iterator v = vectors_begin(); 00841 vectors_iterator v_end = vectors_end(); 00842 unsigned int num = 0; 00843 while((num<vec_num) && (v!=v_end)) 00844 { 00845 num++; 00846 ++v; 00847 } 00848 libmesh_assert (v != v_end); 00849 return *(v->second); 00850 } 00851 00852 00853 00854 const std::string& System::vector_name (const unsigned int vec_num) const 00855 { 00856 const_vectors_iterator v = vectors_begin(); 00857 const_vectors_iterator v_end = vectors_end(); 00858 unsigned int num = 0; 00859 while((num<vec_num) && (v!=v_end)) 00860 { 00861 num++; 00862 ++v; 00863 } 00864 libmesh_assert (v != v_end); 00865 return v->first; 00866 } 00867 00868 const std::string & System::vector_name (const NumericVector<Number> & vec_reference) const 00869 { 00870 const_vectors_iterator v = vectors_begin(); 00871 const_vectors_iterator v_end = vectors_end(); 00872 00873 for(; v != v_end; ++v) 00874 { 00875 // Check if the current vector is the one whose name we want 00876 if(&vec_reference == v->second) 00877 break; // exit loop if it is 00878 } 00879 00880 // Before returning, make sure we didnt loop till the end and not find any match 00881 libmesh_assert (v != v_end); 00882 00883 // Return the string associated with the current vector 00884 return v->first; 00885 } 00886 00887 void System::set_vector_preservation (const std::string &vec_name, 00888 bool preserve) 00889 { 00890 _vector_projections[vec_name] = preserve; 00891 } 00892 00893 00894 00895 bool System::vector_preservation (const std::string &vec_name) const 00896 { 00897 if (_vector_projections.find(vec_name) == _vector_projections.end()) 00898 return false; 00899 00900 return _vector_projections.find(vec_name)->second; 00901 } 00902 00903 00904 00905 NumericVector<Number> & System::add_sensitivity_solution (unsigned int i) 00906 { 00907 std::ostringstream sensitivity_name; 00908 sensitivity_name << "sensitivity_solution" << i; 00909 00910 return this->add_vector(sensitivity_name.str()); 00911 } 00912 00913 00914 00915 NumericVector<Number> & System::get_sensitivity_solution (unsigned int i) 00916 { 00917 std::ostringstream sensitivity_name; 00918 sensitivity_name << "sensitivity_solution" << i; 00919 00920 return this->get_vector(sensitivity_name.str()); 00921 } 00922 00923 00924 00925 const NumericVector<Number> & System::get_sensitivity_solution (unsigned int i) const 00926 { 00927 std::ostringstream sensitivity_name; 00928 sensitivity_name << "sensitivity_solution" << i; 00929 00930 return this->get_vector(sensitivity_name.str()); 00931 } 00932 00933 00934 00935 NumericVector<Number> & System::add_weighted_sensitivity_solution () 00936 { 00937 return this->add_vector("weighted_sensitivity_solution"); 00938 } 00939 00940 00941 00942 NumericVector<Number> & System::get_weighted_sensitivity_solution () 00943 { 00944 return this->get_vector("weighted_sensitivity_solution"); 00945 } 00946 00947 00948 00949 const NumericVector<Number> & System::get_weighted_sensitivity_solution () const 00950 { 00951 return this->get_vector("weighted_sensitivity_solution"); 00952 } 00953 00954 00955 00956 NumericVector<Number> & System::add_adjoint_solution (unsigned int i) 00957 { 00958 std::ostringstream adjoint_name; 00959 adjoint_name << "adjoint_solution" << i; 00960 00961 return this->add_vector(adjoint_name.str()); 00962 } 00963 00964 00965 00966 NumericVector<Number> & System::get_adjoint_solution (unsigned int i) 00967 { 00968 std::ostringstream adjoint_name; 00969 adjoint_name << "adjoint_solution" << i; 00970 00971 return this->get_vector(adjoint_name.str()); 00972 } 00973 00974 00975 00976 const NumericVector<Number> & System::get_adjoint_solution (unsigned int i) const 00977 { 00978 std::ostringstream adjoint_name; 00979 adjoint_name << "adjoint_solution" << i; 00980 00981 return this->get_vector(adjoint_name.str()); 00982 } 00983 00984 00985 00986 NumericVector<Number> & System::add_weighted_sensitivity_adjoint_solution (unsigned int i) 00987 { 00988 std::ostringstream adjoint_name; 00989 adjoint_name << "weighted_sensitivity_adjoint_solution" << i; 00990 00991 return this->add_vector(adjoint_name.str()); 00992 } 00993 00994 00995 00996 NumericVector<Number> & System::get_weighted_sensitivity_adjoint_solution (unsigned int i) 00997 { 00998 std::ostringstream adjoint_name; 00999 adjoint_name << "weighted_sensitivity_adjoint_solution" << i; 01000 01001 return this->get_vector(adjoint_name.str()); 01002 } 01003 01004 01005 01006 const NumericVector<Number> & System::get_weighted_sensitivity_adjoint_solution (unsigned int i) const 01007 { 01008 std::ostringstream adjoint_name; 01009 adjoint_name << "weighted_sensitivity_adjoint_solution" << i; 01010 01011 return this->get_vector(adjoint_name.str()); 01012 } 01013 01014 01015 01016 NumericVector<Number> & System::add_adjoint_rhs (unsigned int i) 01017 { 01018 std::ostringstream adjoint_rhs_name; 01019 adjoint_rhs_name << "adjoint_rhs" << i; 01020 01021 return this->add_vector(adjoint_rhs_name.str(), false); 01022 } 01023 01024 01025 01026 NumericVector<Number> & System::get_adjoint_rhs (unsigned int i) 01027 { 01028 std::ostringstream adjoint_rhs_name; 01029 adjoint_rhs_name << "adjoint_rhs" << i; 01030 01031 return this->get_vector(adjoint_rhs_name.str()); 01032 } 01033 01034 01035 01036 const NumericVector<Number> & System::get_adjoint_rhs (unsigned int i) const 01037 { 01038 std::ostringstream adjoint_rhs_name; 01039 adjoint_rhs_name << "adjoint_rhs" << i; 01040 01041 return this->get_vector(adjoint_rhs_name.str()); 01042 } 01043 01044 01045 01046 NumericVector<Number> & System::add_sensitivity_rhs (unsigned int i) 01047 { 01048 std::ostringstream sensitivity_rhs_name; 01049 sensitivity_rhs_name << "sensitivity_rhs" << i; 01050 01051 return this->add_vector(sensitivity_rhs_name.str(), false); 01052 } 01053 01054 01055 01056 NumericVector<Number> & System::get_sensitivity_rhs (unsigned int i) 01057 { 01058 std::ostringstream sensitivity_rhs_name; 01059 sensitivity_rhs_name << "sensitivity_rhs" << i; 01060 01061 return this->get_vector(sensitivity_rhs_name.str()); 01062 } 01063 01064 01065 01066 const NumericVector<Number> & System::get_sensitivity_rhs (unsigned int i) const 01067 { 01068 std::ostringstream sensitivity_rhs_name; 01069 sensitivity_rhs_name << "sensitivity_rhs" << i; 01070 01071 return this->get_vector(sensitivity_rhs_name.str()); 01072 } 01073 01074 01075 01076 unsigned int System::add_variable (const std::string& var, 01077 const FEType& type, 01078 const std::set<subdomain_id_type> * const active_subdomains) 01079 { 01080 // Make sure the variable isn't there already 01081 // or if it is, that it's the type we want 01082 for (unsigned int v=0; v<this->n_vars(); v++) 01083 if (this->variable_name(v) == var) 01084 { 01085 if (this->variable_type(v) == type) 01086 return _variables[v].number(); 01087 01088 libMesh::err << "ERROR: incompatible variable " 01089 << var 01090 << " has already been added for this system!" 01091 << std::endl; 01092 libmesh_error(); 01093 } 01094 01095 // Optimize for VariableGroups here - if the user is adding multiple 01096 // variables of the same FEType and subdomain restriction, catch 01097 // that here and add them as members of the same VariableGroup. 01098 // 01099 // start by setting this flag to whatever the user has requested 01100 // and then consider the conditions which should negate it. 01101 bool should_be_in_vg = this->identify_variable_groups(); 01102 01103 // No variable groups, nothing to add to 01104 if (!this->n_variable_groups()) 01105 should_be_in_vg = false; 01106 01107 else 01108 { 01109 VariableGroup &vg(_variable_groups.back()); 01110 01111 // get a pointer to their subdomain restriction, if any. 01112 const std::set<subdomain_id_type> * const 01113 their_active_subdomains (vg.implicitly_active() ? 01114 NULL : &vg.active_subdomains()); 01115 01116 // Different types? 01117 if (vg.type() != type) 01118 should_be_in_vg = false; 01119 01120 // they are restricted, we aren't? 01121 if (their_active_subdomains && !active_subdomains) 01122 should_be_in_vg = false; 01123 01124 // they aren't restriced, we are? 01125 if (!their_active_subdomains && active_subdomains) 01126 should_be_in_vg = false; 01127 01128 if (their_active_subdomains && active_subdomains) 01129 // restricted to different sets? 01130 if (*their_active_subdomains != *active_subdomains) 01131 should_be_in_vg = false; 01132 01133 // OK, after all that, append the variable to the vg if none of the conditions 01134 // were violated 01135 if (should_be_in_vg) 01136 { 01137 const unsigned int curr_n_vars = this->n_vars(); 01138 01139 vg.append (var); 01140 01141 _variables.push_back(vg(vg.n_variables()-1)); 01142 _variable_numbers[var] = curr_n_vars; 01143 return curr_n_vars; 01144 } 01145 } 01146 01147 // otherwise, fall back to adding a single variable group 01148 return this->add_variables (std::vector<std::string>(1, var), 01149 type, 01150 active_subdomains); 01151 } 01152 01153 01154 01155 unsigned int System::add_variable (const std::string& var, 01156 const Order order, 01157 const FEFamily family, 01158 const std::set<subdomain_id_type> * const active_subdomains) 01159 { 01160 return this->add_variable(var, 01161 FEType(order, family), 01162 active_subdomains); 01163 } 01164 01165 01166 01167 unsigned int System::add_variables (const std::vector<std::string> &vars, 01168 const FEType& type, 01169 const std::set<subdomain_id_type> * const active_subdomains) 01170 { 01171 // Make sure the variable isn't there already 01172 // or if it is, that it's the type we want 01173 for (unsigned int ov=0; ov<vars.size(); ov++) 01174 for (unsigned int v=0; v<this->n_vars(); v++) 01175 if (this->variable_name(v) == vars[ov]) 01176 { 01177 if (this->variable_type(v) == type) 01178 return _variables[v].number(); 01179 01180 libMesh::err << "ERROR: incompatible variable " 01181 << vars[ov] 01182 << " has already been added for this system!" 01183 << std::endl; 01184 libmesh_error(); 01185 } 01186 01187 const unsigned int curr_n_vars = this->n_vars(); 01188 01189 const unsigned int next_first_component = this->n_components(); 01190 01191 // Add the variable group to the list 01192 _variable_groups.push_back((active_subdomains == NULL) ? 01193 VariableGroup(this, vars, curr_n_vars, 01194 next_first_component, type) : 01195 VariableGroup(this, vars, curr_n_vars, 01196 next_first_component, type, *active_subdomains)); 01197 01198 const VariableGroup &vg (_variable_groups.back()); 01199 01200 // Add each component of the group individually 01201 for (unsigned int v=0; v<vars.size(); v++) 01202 { 01203 _variables.push_back (vg(v)); 01204 _variable_numbers[vars[v]] = curr_n_vars+v; 01205 } 01206 01207 libmesh_assert_equal_to ((curr_n_vars+vars.size()), this->n_vars()); 01208 01209 // BSK - Defer this now to System::init_data() so we can detect 01210 // VariableGroups 12/28/2012 01211 // // Add the variable group to the _dof_map 01212 // _dof_map->add_variable_group (vg); 01213 01214 // Return the number of the new variable 01215 return curr_n_vars+vars.size()-1; 01216 } 01217 01218 01219 01220 unsigned int System::add_variables (const std::vector<std::string> &vars, 01221 const Order order, 01222 const FEFamily family, 01223 const std::set<subdomain_id_type> * const active_subdomains) 01224 { 01225 return this->add_variables(vars, 01226 FEType(order, family), 01227 active_subdomains); 01228 } 01229 01230 01231 01232 bool System::has_variable (const std::string& var) const 01233 { 01234 return _variable_numbers.count(var); 01235 } 01236 01237 01238 01239 unsigned short int System::variable_number (const std::string& var) const 01240 { 01241 // Make sure the variable exists 01242 std::map<std::string, unsigned short int>::const_iterator 01243 pos = _variable_numbers.find(var); 01244 01245 if (pos == _variable_numbers.end()) 01246 { 01247 libMesh::err << "ERROR: variable " 01248 << var 01249 << " does not exist in this system!" 01250 << std::endl; 01251 libmesh_error(); 01252 } 01253 libmesh_assert_equal_to (_variables[pos->second].name(), var); 01254 01255 return pos->second; 01256 } 01257 01258 01259 void System::get_all_variable_numbers(std::vector<unsigned int>& all_variable_numbers) const 01260 { 01261 all_variable_numbers.resize(n_vars()); 01262 01263 // Make sure the variable exists 01264 std::map<std::string, unsigned short int>::const_iterator 01265 it = _variable_numbers.begin(); 01266 std::map<std::string, unsigned short int>::const_iterator 01267 it_end = _variable_numbers.end(); 01268 01269 unsigned int count = 0; 01270 for( ; it != it_end; ++it) 01271 { 01272 all_variable_numbers[count] = it->second; 01273 count++; 01274 } 01275 } 01276 01277 01278 void System::local_dof_indices(const unsigned int var, 01279 std::set<dof_id_type> & var_indices) const 01280 { 01281 // Make sure the set is clear 01282 var_indices.clear(); 01283 01284 std::vector<dof_id_type> dof_indices; 01285 01286 // Begin the loop over the elements 01287 MeshBase::const_element_iterator el = 01288 this->get_mesh().active_local_elements_begin(); 01289 const MeshBase::const_element_iterator end_el = 01290 this->get_mesh().active_local_elements_end(); 01291 01292 const dof_id_type 01293 first_local = this->get_dof_map().first_dof(), 01294 end_local = this->get_dof_map().end_dof(); 01295 01296 for ( ; el != end_el; ++el) 01297 { 01298 const Elem* elem = *el; 01299 this->get_dof_map().dof_indices (elem, dof_indices, var); 01300 01301 for(unsigned int i=0; i<dof_indices.size(); i++) 01302 { 01303 dof_id_type dof = dof_indices[i]; 01304 01305 //If the dof is owned by the local processor 01306 if(first_local <= dof && dof < end_local) 01307 var_indices.insert(dof_indices[i]); 01308 } 01309 } 01310 } 01311 01312 01313 01314 void System::zero_variable (NumericVector<Number>& v, unsigned int var_num) const 01315 { 01316 /* Make sure the call makes sense. */ 01317 libmesh_assert_less (var_num, this->n_vars()); 01318 01319 /* Get a reference to the mesh. */ 01320 const MeshBase& mesh = this->get_mesh(); 01321 01322 /* Check which system we are. */ 01323 const unsigned int sys_num = this->number(); 01324 01325 /* Loop over nodes. */ 01326 { 01327 MeshBase::const_node_iterator it = mesh.local_nodes_begin(); 01328 const MeshBase::const_node_iterator end_it = mesh.local_nodes_end(); 01329 for ( ; it != end_it; ++it) 01330 { 01331 const Node* node = *it; 01332 unsigned int n_comp = node->n_comp(sys_num,var_num); 01333 for(unsigned int i=0; i<n_comp; i++) 01334 { 01335 const dof_id_type index = node->dof_number(sys_num,var_num,i); 01336 v.set(index,0.0); 01337 } 01338 } 01339 } 01340 01341 /* Loop over elements. */ 01342 { 01343 MeshBase::const_element_iterator it = mesh.active_local_elements_begin(); 01344 const MeshBase::const_element_iterator end_it = mesh.active_local_elements_end(); 01345 for ( ; it != end_it; ++it) 01346 { 01347 const Elem* elem = *it; 01348 unsigned int n_comp = elem->n_comp(sys_num,var_num); 01349 for(unsigned int i=0; i<n_comp; i++) 01350 { 01351 const dof_id_type index = elem->dof_number(sys_num,var_num,i); 01352 v.set(index,0.0); 01353 } 01354 } 01355 } 01356 } 01357 01358 01359 01360 Real System::discrete_var_norm(const NumericVector<Number>& v, 01361 unsigned int var, 01362 FEMNormType norm_type) const 01363 { 01364 std::set<dof_id_type> var_indices; 01365 local_dof_indices(var, var_indices); 01366 01367 if(norm_type == DISCRETE_L1) 01368 return v.subset_l1_norm(var_indices); 01369 if(norm_type == DISCRETE_L2) 01370 return v.subset_l2_norm(var_indices); 01371 if(norm_type == DISCRETE_L_INF) 01372 return v.subset_linfty_norm(var_indices); 01373 else 01374 libmesh_error(); 01375 } 01376 01377 01378 01379 Real System::calculate_norm(const NumericVector<Number>& v, 01380 unsigned int var, 01381 FEMNormType norm_type) const 01382 { 01383 //short circuit to save time 01384 if(norm_type == DISCRETE_L1 || 01385 norm_type == DISCRETE_L2 || 01386 norm_type == DISCRETE_L_INF) 01387 return discrete_var_norm(v,var,norm_type); 01388 01389 // Not a discrete norm 01390 std::vector<FEMNormType> norms(this->n_vars(), L2); 01391 std::vector<Real> weights(this->n_vars(), 0.0); 01392 norms[var] = norm_type; 01393 weights[var] = 1.0; 01394 Real val = this->calculate_norm(v, SystemNorm(norms, weights)); 01395 return val; 01396 } 01397 01398 01399 01400 Real System::calculate_norm(const NumericVector<Number>& v, 01401 const SystemNorm &norm) const 01402 { 01403 // This function must be run on all processors at once 01404 parallel_only(); 01405 01406 START_LOG ("calculate_norm()", "System"); 01407 01408 // Zero the norm before summation 01409 Real v_norm = 0.; 01410 01411 if (norm.is_discrete()) 01412 { 01413 STOP_LOG ("calculate_norm()", "System"); 01414 //Check to see if all weights are 1.0 and all types are equal 01415 FEMNormType norm_type0 = norm.type(0); 01416 unsigned int check_var = 0; 01417 for (; check_var != this->n_vars(); ++check_var) 01418 if((norm.weight(check_var) != 1.0) || (norm.type(check_var) != norm_type0)) 01419 break; 01420 01421 //All weights were 1.0 so just do the full vector discrete norm 01422 if(check_var == this->n_vars()) 01423 { 01424 if(norm_type0 == DISCRETE_L1) 01425 return v.l1_norm(); 01426 if(norm_type0 == DISCRETE_L2) 01427 return v.l2_norm(); 01428 if(norm_type0 == DISCRETE_L_INF) 01429 return v.linfty_norm(); 01430 else 01431 libmesh_error(); 01432 } 01433 01434 for (unsigned int var=0; var != this->n_vars(); ++var) 01435 { 01436 // Skip any variables we don't need to integrate 01437 if (norm.weight(var) == 0.0) 01438 continue; 01439 01440 v_norm += norm.weight(var) * discrete_var_norm(v, var, norm.type(var)); 01441 } 01442 01443 return v_norm; 01444 } 01445 01446 // Localize the potentially parallel vector 01447 AutoPtr<NumericVector<Number> > local_v = NumericVector<Number>::build(); 01448 local_v->init(v.size(), true, SERIAL); 01449 v.localize (*local_v, _dof_map->get_send_list()); 01450 01451 unsigned int dim = this->get_mesh().mesh_dimension(); 01452 01453 // I'm not sure how best to mix Hilbert norms on some variables (for 01454 // which we'll want to square then sum then square root) with norms 01455 // like L_inf (for which we'll just want to take an absolute value 01456 // and then sum). 01457 bool using_hilbert_norm = true, 01458 using_nonhilbert_norm = true; 01459 01460 // Loop over all variables 01461 for (unsigned int var=0; var != this->n_vars(); ++var) 01462 { 01463 // Skip any variables we don't need to integrate 01464 Real norm_weight_sq = norm.weight_sq(var); 01465 if (norm_weight_sq == 0.0) 01466 continue; 01467 Real norm_weight = norm.weight(var); 01468 01469 // Check for unimplemented norms (rather than just returning 0). 01470 FEMNormType norm_type = norm.type(var); 01471 if((norm_type==H1) || 01472 (norm_type==H2) || 01473 (norm_type==L2) || 01474 (norm_type==H1_SEMINORM) || 01475 (norm_type==H2_SEMINORM)) 01476 { 01477 if (!using_hilbert_norm) 01478 libmesh_not_implemented(); 01479 using_nonhilbert_norm = false; 01480 } 01481 else if ((norm_type==L1) || 01482 (norm_type==L_INF) || 01483 (norm_type==W1_INF_SEMINORM) || 01484 (norm_type==W2_INF_SEMINORM)) 01485 { 01486 if (!using_nonhilbert_norm) 01487 libmesh_not_implemented(); 01488 using_hilbert_norm = false; 01489 } 01490 else 01491 libmesh_not_implemented(); 01492 01493 const FEType& fe_type = this->get_dof_map().variable_type(var); 01494 AutoPtr<QBase> qrule = 01495 fe_type.default_quadrature_rule (dim); 01496 AutoPtr<FEBase> fe 01497 (FEBase::build(dim, fe_type)); 01498 fe->attach_quadrature_rule (qrule.get()); 01499 01500 const std::vector<Real>& JxW = fe->get_JxW(); 01501 const std::vector<std::vector<Real> >* phi = NULL; 01502 if (norm_type == H1 || 01503 norm_type == H2 || 01504 norm_type == L2 || 01505 norm_type == L1 || 01506 norm_type == L_INF) 01507 phi = &(fe->get_phi()); 01508 01509 const std::vector<std::vector<RealGradient> >* dphi = NULL; 01510 if (norm_type == H1 || 01511 norm_type == H2 || 01512 norm_type == H1_SEMINORM || 01513 norm_type == W1_INF_SEMINORM) 01514 dphi = &(fe->get_dphi()); 01515 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 01516 const std::vector<std::vector<RealTensor> >* d2phi = NULL; 01517 if (norm_type == H2 || 01518 norm_type == H2_SEMINORM || 01519 norm_type == W2_INF_SEMINORM) 01520 d2phi = &(fe->get_d2phi()); 01521 #endif 01522 01523 std::vector<dof_id_type> dof_indices; 01524 01525 // Begin the loop over the elements 01526 MeshBase::const_element_iterator el = 01527 this->get_mesh().active_local_elements_begin(); 01528 const MeshBase::const_element_iterator end_el = 01529 this->get_mesh().active_local_elements_end(); 01530 01531 for ( ; el != end_el; ++el) 01532 { 01533 const Elem* elem = *el; 01534 01535 fe->reinit (elem); 01536 01537 this->get_dof_map().dof_indices (elem, dof_indices, var); 01538 01539 const unsigned int n_qp = qrule->n_points(); 01540 01541 const unsigned int n_sf = libmesh_cast_int<unsigned int> 01542 (dof_indices.size()); 01543 01544 // Begin the loop over the Quadrature points. 01545 for (unsigned int qp=0; qp<n_qp; qp++) 01546 { 01547 if (norm_type == L1) 01548 { 01549 Number u_h = 0.; 01550 for (unsigned int i=0; i != n_sf; ++i) 01551 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]); 01552 v_norm += norm_weight * 01553 JxW[qp] * std::abs(u_h); 01554 } 01555 01556 if (norm_type == L_INF) 01557 { 01558 Number u_h = 0.; 01559 for (unsigned int i=0; i != n_sf; ++i) 01560 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]); 01561 v_norm = std::max(v_norm, norm_weight * std::abs(u_h)); 01562 } 01563 01564 if (norm_type == H1 || 01565 norm_type == H2 || 01566 norm_type == L2) 01567 { 01568 Number u_h = 0.; 01569 for (unsigned int i=0; i != n_sf; ++i) 01570 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]); 01571 v_norm += norm_weight_sq * 01572 JxW[qp] * TensorTools::norm_sq(u_h); 01573 } 01574 01575 if (norm_type == H1 || 01576 norm_type == H2 || 01577 norm_type == H1_SEMINORM) 01578 { 01579 Gradient grad_u_h; 01580 for (unsigned int i=0; i != n_sf; ++i) 01581 grad_u_h.add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i])); 01582 v_norm += norm_weight_sq * 01583 JxW[qp] * grad_u_h.size_sq(); 01584 } 01585 01586 if (norm_type == W1_INF_SEMINORM) 01587 { 01588 Gradient grad_u_h; 01589 for (unsigned int i=0; i != n_sf; ++i) 01590 grad_u_h.add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i])); 01591 v_norm = std::max(v_norm, norm_weight * grad_u_h.size()); 01592 } 01593 01594 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 01595 if (norm_type == H2 || 01596 norm_type == H2_SEMINORM) 01597 { 01598 Tensor hess_u_h; 01599 for (unsigned int i=0; i != n_sf; ++i) 01600 hess_u_h.add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i])); 01601 v_norm += norm_weight_sq * 01602 JxW[qp] * hess_u_h.size_sq(); 01603 } 01604 01605 if (norm_type == W2_INF_SEMINORM) 01606 { 01607 Tensor hess_u_h; 01608 for (unsigned int i=0; i != n_sf; ++i) 01609 hess_u_h.add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i])); 01610 v_norm = std::max(v_norm, norm_weight * hess_u_h.size()); 01611 } 01612 #endif 01613 } 01614 } 01615 } 01616 01617 if (using_hilbert_norm) 01618 { 01619 CommWorld.sum(v_norm); 01620 v_norm = std::sqrt(v_norm); 01621 } 01622 else 01623 { 01624 CommWorld.max(v_norm); 01625 } 01626 01627 STOP_LOG ("calculate_norm()", "System"); 01628 01629 return v_norm; 01630 } 01631 01632 01633 01634 std::string System::get_info() const 01635 { 01636 std::ostringstream oss; 01637 01638 01639 const std::string& sys_name = this->name(); 01640 01641 oss << " System #" << this->number() << ", \"" << sys_name << "\"\n" 01642 << " Type \"" << this->system_type() << "\"\n" 01643 << " Variables="; 01644 01645 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 01646 { 01647 const VariableGroup &vg_description (this->variable_group(vg)); 01648 01649 if (vg_description.n_variables() > 1) oss << "{ "; 01650 for (unsigned int vn=0; vn<vg_description.n_variables(); vn++) 01651 oss << "\"" << vg_description.name(vn) << "\" "; 01652 if (vg_description.n_variables() > 1) oss << "} "; 01653 } 01654 01655 oss << '\n'; 01656 01657 oss << " Finite Element Types="; 01658 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS 01659 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 01660 oss << "\"" 01661 << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().family) 01662 << "\" "; 01663 #else 01664 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 01665 { 01666 oss << "\"" 01667 << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().family) 01668 << "\", \"" 01669 << Utility::enum_to_string<FEFamily>(this->get_dof_map().variable_group(vg).type().radial_family) 01670 << "\" "; 01671 } 01672 01673 oss << '\n' << " Infinite Element Mapping="; 01674 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 01675 oss << "\"" 01676 << Utility::enum_to_string<InfMapType>(this->get_dof_map().variable_group(vg).type().inf_map) 01677 << "\" "; 01678 #endif 01679 01680 oss << '\n'; 01681 01682 oss << " Approximation Orders="; 01683 for (unsigned int vg=0; vg<this->n_variable_groups(); vg++) 01684 { 01685 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS 01686 oss << "\"" 01687 << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().order) 01688 << "\" "; 01689 #else 01690 oss << "\"" 01691 << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().order) 01692 << "\", \"" 01693 << Utility::enum_to_string<Order>(this->get_dof_map().variable_group(vg).type().radial_order) 01694 << "\" "; 01695 #endif 01696 } 01697 01698 oss << '\n'; 01699 01700 oss << " n_dofs()=" << this->n_dofs() << '\n'; 01701 oss << " n_local_dofs()=" << this->n_local_dofs() << '\n'; 01702 #ifdef LIBMESH_ENABLE_CONSTRAINTS 01703 oss << " n_constrained_dofs()=" << this->n_constrained_dofs() << '\n'; 01704 oss << " n_local_constrained_dofs()=" << this->n_local_constrained_dofs() << '\n'; 01705 #endif 01706 01707 oss << " " << "n_vectors()=" << this->n_vectors() << '\n'; 01708 oss << " " << "n_matrices()=" << this->n_matrices() << '\n'; 01709 // oss << " " << "n_additional_matrices()=" << this->n_additional_matrices() << '\n'; 01710 01711 oss << this->get_dof_map().get_info(); 01712 01713 return oss.str(); 01714 } 01715 01716 01717 01718 void System::attach_init_function (void fptr(EquationSystems& es, 01719 const std::string& name)) 01720 { 01721 libmesh_assert(fptr); 01722 01723 if (_init_system_object != NULL) 01724 { 01725 libmesh_here(); 01726 libMesh::out << "WARNING: Cannot specify both initialization function and object!" 01727 << std::endl; 01728 01729 _init_system_object = NULL; 01730 } 01731 01732 _init_system_function = fptr; 01733 } 01734 01735 01736 01737 void System::attach_init_object (System::Initialization& init_in) 01738 { 01739 if (_init_system_function != NULL) 01740 { 01741 libmesh_here(); 01742 libMesh::out << "WARNING: Cannot specify both initialization object and function!" 01743 << std::endl; 01744 01745 _init_system_function = NULL; 01746 } 01747 01748 _init_system_object = &init_in; 01749 } 01750 01751 01752 01753 void System::attach_assemble_function (void fptr(EquationSystems& es, 01754 const std::string& name)) 01755 { 01756 libmesh_assert(fptr); 01757 01758 if (_assemble_system_object != NULL) 01759 { 01760 libmesh_here(); 01761 libMesh::out << "WARNING: Cannot specify both assembly function and object!" 01762 << std::endl; 01763 01764 _assemble_system_object = NULL; 01765 } 01766 01767 _assemble_system_function = fptr; 01768 } 01769 01770 01771 01772 void System::attach_assemble_object (System::Assembly& assemble_in) 01773 { 01774 if (_assemble_system_function != NULL) 01775 { 01776 libmesh_here(); 01777 libMesh::out << "WARNING: Cannot specify both assembly object and function!" 01778 << std::endl; 01779 01780 _assemble_system_function = NULL; 01781 } 01782 01783 _assemble_system_object = &assemble_in; 01784 } 01785 01786 01787 01788 void System::attach_constraint_function(void fptr(EquationSystems& es, 01789 const std::string& name)) 01790 { 01791 libmesh_assert(fptr); 01792 01793 if (_constrain_system_object != NULL) 01794 { 01795 libmesh_here(); 01796 libMesh::out << "WARNING: Cannot specify both constraint function and object!" 01797 << std::endl; 01798 01799 _constrain_system_object = NULL; 01800 } 01801 01802 _constrain_system_function = fptr; 01803 } 01804 01805 01806 01807 void System::attach_constraint_object (System::Constraint& constrain) 01808 { 01809 if (_constrain_system_function != NULL) 01810 { 01811 libmesh_here(); 01812 libMesh::out << "WARNING: Cannot specify both constraint object and function!" 01813 << std::endl; 01814 01815 _constrain_system_function = NULL; 01816 } 01817 01818 _constrain_system_object = &constrain; 01819 } 01820 01821 01822 01823 void System::attach_QOI_function(void fptr(EquationSystems&, 01824 const std::string&, 01825 const QoISet&)) 01826 { 01827 libmesh_assert(fptr); 01828 01829 if (_qoi_evaluate_object != NULL) 01830 { 01831 libmesh_here(); 01832 libMesh::out << "WARNING: Cannot specify both QOI function and object!" 01833 << std::endl; 01834 01835 _qoi_evaluate_object = NULL; 01836 } 01837 01838 _qoi_evaluate_function = fptr; 01839 } 01840 01841 01842 01843 void System::attach_QOI_object (QOI& qoi_in) 01844 { 01845 if (_qoi_evaluate_function != NULL) 01846 { 01847 libmesh_here(); 01848 libMesh::out << "WARNING: Cannot specify both QOI object and function!" 01849 << std::endl; 01850 01851 _qoi_evaluate_function = NULL; 01852 } 01853 01854 _qoi_evaluate_object = &qoi_in; 01855 } 01856 01857 01858 01859 void System::attach_QOI_derivative(void fptr(EquationSystems&, 01860 const std::string&, 01861 const QoISet&)) 01862 { 01863 libmesh_assert(fptr); 01864 01865 if (_qoi_evaluate_derivative_object != NULL) 01866 { 01867 libmesh_here(); 01868 libMesh::out << "WARNING: Cannot specify both QOI derivative function and object!" 01869 << std::endl; 01870 01871 _qoi_evaluate_derivative_object = NULL; 01872 } 01873 01874 _qoi_evaluate_derivative_function = fptr; 01875 } 01876 01877 01878 01879 void System::attach_QOI_derivative_object (QOIDerivative& qoi_derivative) 01880 { 01881 if (_qoi_evaluate_derivative_function != NULL) 01882 { 01883 libmesh_here(); 01884 libMesh::out << "WARNING: Cannot specify both QOI derivative object and function!" 01885 << std::endl; 01886 01887 _qoi_evaluate_derivative_function = NULL; 01888 } 01889 01890 _qoi_evaluate_derivative_object = &qoi_derivative; 01891 } 01892 01893 01894 01895 void System::user_initialization () 01896 { 01897 // Call the user-provided intialization function, 01898 // if it was provided 01899 if (_init_system_function != NULL) 01900 this->_init_system_function (_equation_systems, this->name()); 01901 01902 // ...or the user-provided initialization object. 01903 else if (_init_system_object != NULL) 01904 this->_init_system_object->initialize(); 01905 } 01906 01907 01908 01909 void System::user_assembly () 01910 { 01911 // Call the user-provided assembly function, 01912 // if it was provided 01913 if (_assemble_system_function != NULL) 01914 this->_assemble_system_function (_equation_systems, this->name()); 01915 01916 // ...or the user-provided assembly object. 01917 else if (_assemble_system_object != NULL) 01918 this->_assemble_system_object->assemble(); 01919 } 01920 01921 01922 01923 void System::user_constrain () 01924 { 01925 // Call the user-provided constraint function, 01926 // if it was provided 01927 if (_constrain_system_function!= NULL) 01928 this->_constrain_system_function(_equation_systems, this->name()); 01929 01930 // ...or the user-provided constraint object. 01931 else if (_constrain_system_object != NULL) 01932 this->_constrain_system_object->constrain(); 01933 } 01934 01935 01936 01937 void System::user_QOI (const QoISet& qoi_indices) 01938 { 01939 // Call the user-provided quantity of interest function, 01940 // if it was provided 01941 if (_qoi_evaluate_function != NULL) 01942 this->_qoi_evaluate_function(_equation_systems, this->name(), qoi_indices); 01943 01944 // ...or the user-provided QOI function object. 01945 else if (_qoi_evaluate_object != NULL) 01946 this->_qoi_evaluate_object->qoi(qoi_indices); 01947 } 01948 01949 01950 01951 void System::user_QOI_derivative (const QoISet& qoi_indices) 01952 { 01953 // Call the user-provided quantity of interest derivative, 01954 // if it was provided 01955 if (_qoi_evaluate_derivative_function != NULL) 01956 this->_qoi_evaluate_derivative_function(_equation_systems, this->name(), qoi_indices); 01957 01958 // ...or the user-provided QOI derivative function object. 01959 else if (_qoi_evaluate_derivative_object != NULL) 01960 this->_qoi_evaluate_derivative_object->qoi_derivative(qoi_indices); 01961 } 01962 01963 01964 01965 Number System::point_value(unsigned int var, const Point &p, const bool insist_on_success) const 01966 { 01967 // This function must be called on every processor; there's no 01968 // telling where in the partition p falls. 01969 parallel_only(); 01970 01971 // And every processor had better agree about which point we're 01972 // looking for 01973 #ifndef NDEBUG 01974 CommWorld.verify(p); 01975 #endif // NDEBUG 01976 01977 // Get a reference to the mesh object associated with the system object that calls this function 01978 const MeshBase &mesh = this->get_mesh(); 01979 01980 // Use an existing PointLocator or create a new one 01981 AutoPtr<PointLocatorBase> locator_ptr = mesh.sub_point_locator(); 01982 PointLocatorBase& locator = *locator_ptr; 01983 01984 if (!insist_on_success) 01985 locator.enable_out_of_mesh_mode(); 01986 01987 // Get a pointer to the element that contains P 01988 const Elem *e = locator(p); 01989 01990 Number u = 0; 01991 01992 if (e && e->processor_id() == libMesh::processor_id()) 01993 u = point_value(var, p, *e); 01994 01995 // If I have an element containing p, then let's let everyone know 01996 processor_id_type lowest_owner = 01997 (e && (e->processor_id() == libMesh::processor_id())) ? 01998 libMesh::processor_id() : libMesh::n_processors(); 01999 CommWorld.min(lowest_owner); 02000 02001 // Everybody should get their value from a processor that was able 02002 // to compute it. 02003 // If nobody admits owning the point, we have a problem. 02004 if (lowest_owner != libMesh::n_processors()) 02005 CommWorld.broadcast(u, lowest_owner); 02006 else 02007 libmesh_assert(!insist_on_success); 02008 02009 return u; 02010 } 02011 02012 Number System::point_value(unsigned int var, const Point &p, const Elem &e) const 02013 { 02014 libmesh_assert_equal_to (e.processor_id(), libMesh::processor_id()); 02015 02016 // Ensuring that the given point is really in the element is an 02017 // expensive assert, but as long as debugging is turned on we might 02018 // as well try to catch a particularly nasty potential error 02019 libmesh_assert (e.contains_point(p)); 02020 02021 // Get the dof map to get the proper indices for our computation 02022 const DofMap& dof_map = this->get_dof_map(); 02023 02024 // Need dof_indices for phi[i][j] 02025 std::vector<dof_id_type> dof_indices; 02026 02027 // Fill in the dof_indices for our element 02028 dof_map.dof_indices (&e, dof_indices, var); 02029 02030 // Get the no of dofs assciated with this point 02031 const unsigned int num_dofs = libmesh_cast_int<unsigned int> 02032 (dof_indices.size()); 02033 02034 FEType fe_type = dof_map.variable_type(0); 02035 02036 // Build a FE so we can calculate u(p) 02037 AutoPtr<FEBase> fe (FEBase::build(e.dim(), fe_type)); 02038 02039 // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h 02040 // Build a vector of point co-ordinates to send to reinit 02041 std::vector<Point> coor(1, FEInterface::inverse_map(e.dim(), fe_type, &e, p)); 02042 02043 // Get the shape function values 02044 const std::vector<std::vector<Real> >& phi = fe->get_phi(); 02045 02046 // Reinitialize the element and compute the shape function values at coor 02047 fe->reinit (&e, &coor); 02048 02049 // Get ready to accumulate a value 02050 Number u = 0; 02051 02052 for (unsigned int l=0; l<num_dofs; l++) 02053 { 02054 u += phi[l][0]*this->current_solution (dof_indices[l]); 02055 } 02056 02057 return u; 02058 } 02059 02060 02061 02062 Gradient System::point_gradient(unsigned int var, const Point &p, const bool insist_on_success) const 02063 { 02064 // This function must be called on every processor; there's no 02065 // telling where in the partition p falls. 02066 parallel_only(); 02067 02068 // And every processor had better agree about which point we're 02069 // looking for 02070 #ifndef NDEBUG 02071 CommWorld.verify(p); 02072 #endif // NDEBUG 02073 02074 // Get a reference to the mesh object associated with the system object that calls this function 02075 const MeshBase &mesh = this->get_mesh(); 02076 02077 // Use an existing PointLocator or create a new one 02078 AutoPtr<PointLocatorBase> locator_ptr = mesh.sub_point_locator(); 02079 PointLocatorBase& locator = *locator_ptr; 02080 02081 if (!insist_on_success) 02082 locator.enable_out_of_mesh_mode(); 02083 02084 // Get a pointer to the element that contains P 02085 const Elem *e = locator(p); 02086 02087 Gradient grad_u; 02088 02089 if (e && e->processor_id() == libMesh::processor_id()) 02090 grad_u = point_gradient(var, p, *e); 02091 02092 // If I have an element containing p, then let's let everyone know 02093 processor_id_type lowest_owner = 02094 (e && (e->processor_id() == libMesh::processor_id())) ? 02095 libMesh::processor_id() : libMesh::n_processors(); 02096 CommWorld.min(lowest_owner); 02097 02098 // Everybody should get their value from a processor that was able 02099 // to compute it. 02100 // If nobody admits owning the point, we may have a problem. 02101 if (lowest_owner != libMesh::n_processors()) 02102 CommWorld.broadcast(grad_u, lowest_owner); 02103 else 02104 libmesh_assert(!insist_on_success); 02105 02106 return grad_u; 02107 } 02108 02109 02110 Gradient System::point_gradient(unsigned int var, const Point &p, const Elem &e) const 02111 { 02112 libmesh_assert_equal_to (e.processor_id(), libMesh::processor_id()); 02113 02114 // Ensuring that the given point is really in the element is an 02115 // expensive assert, but as long as debugging is turned on we might 02116 // as well try to catch a particularly nasty potential error 02117 libmesh_assert (e.contains_point(p)); 02118 02119 // Get the dof map to get the proper indices for our computation 02120 const DofMap& dof_map = this->get_dof_map(); 02121 02122 // Need dof_indices for phi[i][j] 02123 std::vector<dof_id_type> dof_indices; 02124 02125 // Fill in the dof_indices for our element 02126 dof_map.dof_indices (&e, dof_indices, var); 02127 02128 // Get the no of dofs assciated with this point 02129 const unsigned int num_dofs = libmesh_cast_int<unsigned int> 02130 (dof_indices.size()); 02131 02132 FEType fe_type = dof_map.variable_type(0); 02133 02134 // Build a FE again so we can calculate u(p) 02135 AutoPtr<FEBase> fe (FEBase::build(e.dim(), fe_type)); 02136 02137 // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h 02138 // Build a vector of point co-ordinates to send to reinit 02139 std::vector<Point> coor(1, FEInterface::inverse_map(e.dim(), fe_type, &e, p)); 02140 02141 // Get the values of the shape function derivatives 02142 const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi(); 02143 02144 // Reinitialize the element and compute the shape function values at coor 02145 fe->reinit (&e, &coor); 02146 02147 // Get ready to accumulate a gradient 02148 Gradient grad_u; 02149 02150 for (unsigned int l=0; l<num_dofs; l++) 02151 { 02152 grad_u.add_scaled (dphi[l][0], this->current_solution (dof_indices[l])); 02153 } 02154 02155 return grad_u; 02156 } 02157 02158 02159 // We can only accumulate a hessian with --enable-second 02160 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 02161 Tensor System::point_hessian(unsigned int var, const Point &p, const bool insist_on_success) const 02162 { 02163 // This function must be called on every processor; there's no 02164 // telling where in the partition p falls. 02165 parallel_only(); 02166 02167 // And every processor had better agree about which point we're 02168 // looking for 02169 #ifndef NDEBUG 02170 CommWorld.verify(p); 02171 #endif // NDEBUG 02172 02173 // Get a reference to the mesh object associated with the system object that calls this function 02174 const MeshBase &mesh = this->get_mesh(); 02175 02176 // Use an existing PointLocator or create a new one 02177 AutoPtr<PointLocatorBase> locator_ptr = mesh.sub_point_locator(); 02178 PointLocatorBase& locator = *locator_ptr; 02179 02180 if (!insist_on_success) 02181 locator.enable_out_of_mesh_mode(); 02182 02183 // Get a pointer to the element that contains P 02184 const Elem *e = locator(p); 02185 02186 Tensor hess_u; 02187 02188 if (e && e->processor_id() == libMesh::processor_id()) 02189 hess_u = point_hessian(var, p, *e); 02190 02191 // If I have an element containing p, then let's let everyone know 02192 processor_id_type lowest_owner = 02193 (e && (e->processor_id() == libMesh::processor_id())) ? 02194 libMesh::processor_id() : libMesh::n_processors(); 02195 CommWorld.min(lowest_owner); 02196 02197 // Everybody should get their value from a processor that was able 02198 // to compute it. 02199 // If nobody admits owning the point, we may have a problem. 02200 if (lowest_owner != libMesh::n_processors()) 02201 CommWorld.broadcast(hess_u, lowest_owner); 02202 else 02203 libmesh_assert(!insist_on_success); 02204 02205 return hess_u; 02206 } 02207 02208 Tensor System::point_hessian(unsigned int var, const Point &p, const Elem &e) const 02209 { 02210 libmesh_assert_equal_to (e.processor_id(), libMesh::processor_id()); 02211 02212 // Ensuring that the given point is really in the element is an 02213 // expensive assert, but as long as debugging is turned on we might 02214 // as well try to catch a particularly nasty potential error 02215 libmesh_assert (e.contains_point(p)); 02216 02217 // Get the dof map to get the proper indices for our computation 02218 const DofMap& dof_map = this->get_dof_map(); 02219 02220 // Need dof_indices for phi[i][j] 02221 std::vector<dof_id_type> dof_indices; 02222 02223 // Fill in the dof_indices for our element 02224 dof_map.dof_indices (&e, dof_indices, var); 02225 02226 // Get the no of dofs assciated with this point 02227 const unsigned int num_dofs = libmesh_cast_int<unsigned int> 02228 (dof_indices.size()); 02229 02230 FEType fe_type = dof_map.variable_type(0); 02231 02232 // Build a FE again so we can calculate u(p) 02233 AutoPtr<FEBase> fe (FEBase::build(e.dim(), fe_type)); 02234 02235 // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h 02236 // Build a vector of point co-ordinates to send to reinit 02237 std::vector<Point> coor(1, FEInterface::inverse_map(e.dim(), fe_type, &e, p)); 02238 02239 // Get the values of the shape function derivatives 02240 const std::vector<std::vector<RealTensor> >& d2phi = fe->get_d2phi(); 02241 02242 // Reinitialize the element and compute the shape function values at coor 02243 fe->reinit (&e, &coor); 02244 02245 // Get ready to accumulate a hessian 02246 Tensor hess_u; 02247 02248 for (unsigned int l=0; l<num_dofs; l++) 02249 { 02250 hess_u.add_scaled (d2phi[l][0], this->current_solution (dof_indices[l])); 02251 } 02252 02253 return hess_u; 02254 } 02255 #else 02256 Tensor System::point_hessian(unsigned int, const Point &, const bool) const 02257 { 02258 // We can only accumulate a hessian with --enable-second 02259 libmesh_error(); 02260 02261 // Avoid compiler warnings 02262 return Tensor(); 02263 } 02264 02265 Tensor System::point_hessian(unsigned int var, const Point &p, const Elem &e) const 02266 { 02267 // We can only accumulate a hessian with --enable-second 02268 libmesh_error(); 02269 02270 // Avoid compiler warnings 02271 return Tensor(); 02272 } 02273 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 02274 02275 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: