system_projection.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 <vector> 00022 00023 // Local includes 00024 #include "libmesh/boundary_info.h" 00025 #include "libmesh/dense_matrix.h" 00026 #include "libmesh/dense_vector.h" 00027 #include "libmesh/dirichlet_boundaries.h" 00028 #include "libmesh/dof_map.h" 00029 #include "libmesh/elem.h" 00030 #include "libmesh/equation_systems.h" 00031 #include "libmesh/fe_base.h" 00032 #include "libmesh/fe_interface.h" 00033 #include "libmesh/libmesh_logging.h" 00034 #include "libmesh/mesh_base.h" 00035 #include "libmesh/numeric_vector.h" 00036 #include "libmesh/quadrature_gauss.h" 00037 #include "libmesh/system.h" 00038 #include "libmesh/threads.h" 00039 #include "libmesh/wrapped_function.h" 00040 00041 namespace libMesh 00042 { 00043 00044 // ------------------------------------------------------------ 00045 // Helper class definitions 00046 00052 class ProjectVector 00053 { 00054 private: 00055 const System &system; 00056 const NumericVector<Number> &old_vector; 00057 NumericVector<Number> &new_vector; 00058 00059 public: 00060 ProjectVector (const System &system_in, 00061 const NumericVector<Number> &old_v_in, 00062 NumericVector<Number> &new_v_in) : 00063 system(system_in), 00064 old_vector(old_v_in), 00065 new_vector(new_v_in) 00066 {} 00067 00068 void operator()(const ConstElemRange &range) const; 00069 }; 00070 00071 00081 class BuildProjectionList 00082 { 00083 private: 00084 const System &system; 00085 00086 public: 00087 BuildProjectionList (const System &system_in) : 00088 system(system_in), 00089 send_list() 00090 {} 00091 00092 BuildProjectionList (BuildProjectionList &other, Threads::split) : 00093 system(other.system), 00094 send_list() 00095 {} 00096 00097 void unique(); 00098 void operator()(const ConstElemRange &range); 00099 void join (const BuildProjectionList &other); 00100 std::vector<dof_id_type> send_list; 00101 }; 00102 00103 00104 00110 class ProjectSolution 00111 { 00112 private: 00113 const System &system; 00114 00115 AutoPtr<FunctionBase<Number> > f; 00116 AutoPtr<FunctionBase<Gradient> > g; 00117 const Parameters ¶meters; 00118 NumericVector<Number> &new_vector; 00119 00120 public: 00121 ProjectSolution (const System &system_in, 00122 FunctionBase<Number>* f_in, 00123 FunctionBase<Gradient>* g_in, 00124 const Parameters ¶meters_in, 00125 NumericVector<Number> &new_v_in) : 00126 system(system_in), 00127 f(f_in ? f_in->clone() : AutoPtr<FunctionBase<Number> >(NULL)), 00128 g(g_in ? g_in->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)), 00129 parameters(parameters_in), 00130 new_vector(new_v_in) 00131 { 00132 libmesh_assert(f.get()); 00133 f->init(); 00134 if (g.get()) 00135 g->init(); 00136 } 00137 00138 ProjectSolution (const ProjectSolution &in) : 00139 system(in.system), 00140 f(in.f.get() ? in.f->clone() : AutoPtr<FunctionBase<Number> >(NULL)), 00141 g(in.g.get() ? in.g->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)), 00142 parameters(in.parameters), 00143 new_vector(in.new_vector) 00144 { 00145 libmesh_assert(f.get()); 00146 f->init(); 00147 if (g.get()) 00148 g->init(); 00149 } 00150 00151 void operator()(const ConstElemRange &range) const; 00152 }; 00153 00154 00160 class ProjectFEMSolution 00161 { 00162 private: 00163 const System &system; 00164 00165 AutoPtr<FEMFunctionBase<Number> > f; 00166 AutoPtr<FEMFunctionBase<Gradient> > g; 00167 NumericVector<Number> &new_vector; 00168 00169 public: 00170 ProjectFEMSolution (const System &system_in, 00171 FEMFunctionBase<Number>* f_in, 00172 FEMFunctionBase<Gradient>* g_in, 00173 NumericVector<Number> &new_v_in) : 00174 system(system_in), 00175 f(f_in ? f_in->clone() : AutoPtr<FEMFunctionBase<Number> >(NULL)), 00176 g(g_in ? g_in->clone() : AutoPtr<FEMFunctionBase<Gradient> >(NULL)), 00177 new_vector(new_v_in) 00178 { 00179 libmesh_assert(f.get()); 00180 } 00181 00182 ProjectFEMSolution (const ProjectFEMSolution &in) : 00183 system(in.system), 00184 f(in.f.get() ? in.f->clone() : AutoPtr<FEMFunctionBase<Number> >(NULL)), 00185 g(in.g.get() ? in.g->clone() : AutoPtr<FEMFunctionBase<Gradient> >(NULL)), 00186 new_vector(in.new_vector) 00187 { 00188 libmesh_assert(f.get()); 00189 } 00190 00191 void operator()(const ConstElemRange &range) const; 00192 }; 00193 00194 00200 class BoundaryProjectSolution 00201 { 00202 private: 00203 const std::set<boundary_id_type> &b; 00204 const std::vector<unsigned int> &variables; 00205 const System &system; 00206 AutoPtr<FunctionBase<Number> > f; 00207 AutoPtr<FunctionBase<Gradient> > g; 00208 const Parameters ¶meters; 00209 NumericVector<Number> &new_vector; 00210 00211 public: 00212 BoundaryProjectSolution (const std::set<boundary_id_type> &b_in, 00213 const std::vector<unsigned int> &variables_in, 00214 const System &system_in, 00215 FunctionBase<Number>* f_in, 00216 FunctionBase<Gradient>* g_in, 00217 const Parameters ¶meters_in, 00218 NumericVector<Number> &new_v_in) : 00219 b(b_in), 00220 variables(variables_in), 00221 system(system_in), 00222 f(f_in ? f_in->clone() : AutoPtr<FunctionBase<Number> >(NULL)), 00223 g(g_in ? g_in->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)), 00224 parameters(parameters_in), 00225 new_vector(new_v_in) 00226 { 00227 libmesh_assert(f.get()); 00228 f->init(); 00229 if (g.get()) 00230 g->init(); 00231 } 00232 00233 BoundaryProjectSolution (const BoundaryProjectSolution &in) : 00234 b(in.b), 00235 variables(in.variables), 00236 system(in.system), 00237 f(in.f.get() ? in.f->clone() : AutoPtr<FunctionBase<Number> >(NULL)), 00238 g(in.g.get() ? in.g->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)), 00239 parameters(in.parameters), 00240 new_vector(in.new_vector) 00241 { 00242 libmesh_assert(f.get()); 00243 f->init(); 00244 if (g.get()) 00245 g->init(); 00246 } 00247 00248 void operator()(const ConstElemRange &range) const; 00249 }; 00250 00251 00252 00253 // ------------------------------------------------------------ 00254 // System implementation 00255 void System::project_vector (NumericVector<Number>& vector) const 00256 { 00257 // Create a copy of the vector, which currently 00258 // contains the old data. 00259 AutoPtr<NumericVector<Number> > 00260 old_vector (vector.clone()); 00261 00262 // Project the old vector to the new vector 00263 this->project_vector (*old_vector, vector); 00264 } 00265 00266 00272 void System::project_vector (const NumericVector<Number>& old_v, 00273 NumericVector<Number>& new_v) const 00274 { 00275 START_LOG ("project_vector()", "System"); 00276 00283 new_v.clear(); 00284 00285 #ifdef LIBMESH_ENABLE_AMR 00286 00287 // Resize the new vector and get a serial version. 00288 NumericVector<Number> *new_vector_ptr = NULL; 00289 AutoPtr<NumericVector<Number> > new_vector_built; 00290 NumericVector<Number> *local_old_vector; 00291 AutoPtr<NumericVector<Number> > local_old_vector_built; 00292 const NumericVector<Number> *old_vector_ptr = NULL; 00293 00294 ConstElemRange active_local_elem_range 00295 (this->get_mesh().active_local_elements_begin(), 00296 this->get_mesh().active_local_elements_end()); 00297 00298 // If the old vector was uniprocessor, make the new 00299 // vector uniprocessor 00300 if (old_v.type() == SERIAL) 00301 { 00302 new_v.init (this->n_dofs(), false, SERIAL); 00303 new_vector_ptr = &new_v; 00304 old_vector_ptr = &old_v; 00305 } 00306 00307 // Otherwise it is a parallel, distributed vector, which 00308 // we need to localize. 00309 else if (old_v.type() == PARALLEL) 00310 { 00311 // Build a send list for efficient localization 00312 BuildProjectionList projection_list(*this); 00313 Threads::parallel_reduce (active_local_elem_range, 00314 projection_list); 00315 00316 // Create a sorted, unique send_list 00317 projection_list.unique(); 00318 00319 new_v.init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL); 00320 new_vector_built = NumericVector<Number>::build(); 00321 local_old_vector_built = NumericVector<Number>::build(); 00322 new_vector_ptr = new_vector_built.get(); 00323 local_old_vector = local_old_vector_built.get(); 00324 new_vector_ptr->init(this->n_dofs(), false, SERIAL); 00325 local_old_vector->init(old_v.size(), false, SERIAL); 00326 old_v.localize(*local_old_vector, projection_list.send_list); 00327 local_old_vector->close(); 00328 old_vector_ptr = local_old_vector; 00329 } 00330 else if (old_v.type() == GHOSTED) 00331 { 00332 // Build a send list for efficient localization 00333 BuildProjectionList projection_list(*this); 00334 Threads::parallel_reduce (active_local_elem_range, 00335 projection_list); 00336 00337 // Create a sorted, unique send_list 00338 projection_list.unique(); 00339 00340 new_v.init (this->n_dofs(), this->n_local_dofs(), 00341 this->get_dof_map().get_send_list(), false, GHOSTED); 00342 00343 local_old_vector_built = NumericVector<Number>::build(); 00344 new_vector_ptr = &new_v; 00345 local_old_vector = local_old_vector_built.get(); 00346 local_old_vector->init(old_v.size(), old_v.local_size(), 00347 projection_list.send_list, false, GHOSTED); 00348 old_v.localize(*local_old_vector, projection_list.send_list); 00349 local_old_vector->close(); 00350 old_vector_ptr = local_old_vector; 00351 } 00352 else // unknown old_v.type() 00353 { 00354 libMesh::err << "ERROR: Unknown old_v.type() == " << old_v.type() 00355 << std::endl; 00356 libmesh_error(); 00357 } 00358 00359 // Note that the above will have zeroed the new_vector. 00360 // Just to be sure, assert that new_vector_ptr and old_vector_ptr 00361 // were successfully set before trying to deref them. 00362 libmesh_assert(new_vector_ptr); 00363 libmesh_assert(old_vector_ptr); 00364 00365 NumericVector<Number> &new_vector = *new_vector_ptr; 00366 const NumericVector<Number> &old_vector = *old_vector_ptr; 00367 00368 Threads::parallel_for (active_local_elem_range, 00369 ProjectVector(*this, 00370 old_vector, 00371 new_vector) 00372 ); 00373 00374 // Copy the SCALAR dofs from old_vector to new_vector 00375 // Note: We assume that all SCALAR dofs are on the 00376 // processor with highest ID 00377 if(libMesh::processor_id() == (libMesh::n_processors()-1)) 00378 { 00379 const DofMap& dof_map = this->get_dof_map(); 00380 for (unsigned int var=0; var<this->n_vars(); var++) 00381 if(this->variable(var).type().family == SCALAR) 00382 { 00383 // We can just map SCALAR dofs directly across 00384 std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices; 00385 dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false); 00386 dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true); 00387 const unsigned int new_n_dofs = 00388 libmesh_cast_int<unsigned int>(new_SCALAR_indices.size()); 00389 00390 for (unsigned int i=0; i<new_n_dofs; i++) 00391 { 00392 new_vector.set( new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]) ); 00393 } 00394 } 00395 } 00396 00397 new_vector.close(); 00398 00399 // If the old vector was serial, we probably need to send our values 00400 // to other processors 00401 // 00402 // FIXME: I'm not sure how to make a NumericVector do that without 00403 // creating a temporary parallel vector to use localize! - RHS 00404 if (old_v.type() == SERIAL) 00405 { 00406 AutoPtr<NumericVector<Number> > dist_v = NumericVector<Number>::build(); 00407 dist_v->init(this->n_dofs(), this->n_local_dofs(), false, PARALLEL); 00408 dist_v->close(); 00409 00410 for (dof_id_type i=0; i!=dist_v->size(); i++) 00411 if (new_vector(i) != 0.0) 00412 dist_v->set(i, new_vector(i)); 00413 00414 dist_v->close(); 00415 00416 dist_v->localize (new_v, this->get_dof_map().get_send_list()); 00417 new_v.close(); 00418 } 00419 // If the old vector was parallel, we need to update it 00420 // and free the localized copies 00421 else if (old_v.type() == PARALLEL) 00422 { 00423 // We may have to set dof values that this processor doesn't 00424 // own in certain special cases, like LAGRANGE FIRST or 00425 // HERMITE THIRD elements on second-order meshes 00426 for (dof_id_type i=0; i!=new_v.size(); i++) 00427 if (new_vector(i) != 0.0) 00428 new_v.set(i, new_vector(i)); 00429 new_v.close(); 00430 } 00431 00432 this->get_dof_map().enforce_constraints_exactly(*this, &new_v); 00433 00434 #else 00435 00436 // AMR is disabled: simply copy the vector 00437 new_v = old_v; 00438 00439 #endif // #ifdef LIBMESH_ENABLE_AMR 00440 00441 STOP_LOG("project_vector()", "System"); 00442 } 00443 00444 00445 00450 void System::project_solution (Number fptr(const Point& p, 00451 const Parameters& parameters, 00452 const std::string& sys_name, 00453 const std::string& unknown_name), 00454 Gradient gptr(const Point& p, 00455 const Parameters& parameters, 00456 const std::string& sys_name, 00457 const std::string& unknown_name), 00458 const Parameters& parameters) const 00459 { 00460 WrappedFunction<Number> f(*this, fptr, ¶meters); 00461 WrappedFunction<Gradient> g(*this, gptr, ¶meters); 00462 this->project_solution(&f, &g); 00463 } 00464 00465 00470 void System::project_solution (FunctionBase<Number> *f, 00471 FunctionBase<Gradient> *g) const 00472 { 00473 this->project_vector(*solution, f, g); 00474 00475 solution->localize(*current_local_solution, _dof_map->get_send_list()); 00476 } 00477 00478 00483 void System::project_solution (FEMFunctionBase<Number> *f, 00484 FEMFunctionBase<Gradient> *g) const 00485 { 00486 this->project_vector(*solution, f, g); 00487 00488 solution->localize(*current_local_solution, _dof_map->get_send_list()); 00489 } 00490 00491 00496 void System::project_vector (Number fptr(const Point& p, 00497 const Parameters& parameters, 00498 const std::string& sys_name, 00499 const std::string& unknown_name), 00500 Gradient gptr(const Point& p, 00501 const Parameters& parameters, 00502 const std::string& sys_name, 00503 const std::string& unknown_name), 00504 const Parameters& parameters, 00505 NumericVector<Number>& new_vector) const 00506 { 00507 WrappedFunction<Number> f(*this, fptr, ¶meters); 00508 WrappedFunction<Gradient> g(*this, gptr, ¶meters); 00509 this->project_vector(new_vector, &f, &g); 00510 } 00511 00516 void System::project_vector (NumericVector<Number>& new_vector, 00517 FunctionBase<Number> *f, 00518 FunctionBase<Gradient> *g) const 00519 { 00520 START_LOG ("project_vector()", "System"); 00521 00522 Threads::parallel_for 00523 (ConstElemRange (this->get_mesh().active_local_elements_begin(), 00524 this->get_mesh().active_local_elements_end() ), 00525 ProjectSolution(*this, f, g, 00526 this->get_equation_systems().parameters, 00527 new_vector) 00528 ); 00529 00530 // Also, load values into the SCALAR dofs 00531 // Note: We assume that all SCALAR dofs are on the 00532 // processor with highest ID 00533 if(libMesh::processor_id() == (libMesh::n_processors()-1)) 00534 { 00535 // We get different scalars as different 00536 // components from a new-style f functor. 00537 DenseVector<Number> fout(this->n_components()); 00538 bool filled_fout = false; 00539 00540 const DofMap& dof_map = this->get_dof_map(); 00541 for (unsigned int var=0; var<this->n_vars(); var++) 00542 if(this->variable(var).type().family == SCALAR) 00543 { 00544 if (!filled_fout) 00545 { 00546 (*f) (Point(), this->time, fout); 00547 filled_fout = true; 00548 } 00549 00550 std::vector<dof_id_type> SCALAR_indices; 00551 dof_map.SCALAR_dof_indices (SCALAR_indices, var); 00552 const unsigned int n_SCALAR_dofs = 00553 libmesh_cast_int<unsigned int>(SCALAR_indices.size()); 00554 00555 for (unsigned int i=0; i<n_SCALAR_dofs; i++) 00556 { 00557 const dof_id_type global_index = SCALAR_indices[i]; 00558 const unsigned int component_index = 00559 this->variable_scalar_number(var,i); 00560 new_vector.set(global_index, fout(component_index)); 00561 } 00562 } 00563 } 00564 00565 new_vector.close(); 00566 00567 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00568 this->get_dof_map().enforce_constraints_exactly(*this, &new_vector); 00569 #endif 00570 00571 STOP_LOG("project_vector()", "System"); 00572 } 00573 00574 00579 void System::project_vector (NumericVector<Number>& new_vector, 00580 FEMFunctionBase<Number> *f, 00581 FEMFunctionBase<Gradient> *g) const 00582 { 00583 START_LOG ("project_fem_vector()", "System"); 00584 00585 Threads::parallel_for 00586 (ConstElemRange (this->get_mesh().active_local_elements_begin(), 00587 this->get_mesh().active_local_elements_end() ), 00588 ProjectFEMSolution(*this, f, g, new_vector) 00589 ); 00590 00591 // Also, load values into the SCALAR dofs 00592 // Note: We assume that all SCALAR dofs are on the 00593 // processor with highest ID 00594 if(libMesh::processor_id() == (libMesh::n_processors()-1)) 00595 { 00596 // FIXME: Do we want to first check for SCALAR vars before building this? [PB] 00597 FEMContext context( *this ); 00598 00599 const DofMap& dof_map = this->get_dof_map(); 00600 for (unsigned int var=0; var<this->n_vars(); var++) 00601 if(this->variable(var).type().family == SCALAR) 00602 { 00603 // FIXME: We reinit with an arbitrary element in case the user 00604 // doesn't override FEMFunctionBase::component. Is there 00605 // any use case we're missing? [PB] 00606 Elem *el = const_cast<Elem *>(*(this->get_mesh().active_local_elements_begin())); 00607 context.pre_fe_reinit( *this, el ); 00608 00609 std::vector<dof_id_type> SCALAR_indices; 00610 dof_map.SCALAR_dof_indices (SCALAR_indices, var); 00611 const unsigned int n_SCALAR_dofs = 00612 libmesh_cast_int<unsigned int>(SCALAR_indices.size()); 00613 00614 for (unsigned int i=0; i<n_SCALAR_dofs; i++) 00615 { 00616 const dof_id_type global_index = SCALAR_indices[i]; 00617 const unsigned int component_index = 00618 this->variable_scalar_number(var,i); 00619 00620 new_vector.set(global_index, f->component(context, component_index, Point(), this->time)); 00621 } 00622 } 00623 } 00624 00625 new_vector.close(); 00626 00627 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00628 this->get_dof_map().enforce_constraints_exactly(*this, &new_vector); 00629 #endif 00630 00631 STOP_LOG("project_fem_vector()", "System"); 00632 } 00633 00634 00640 void System::boundary_project_solution 00641 (const std::set<boundary_id_type> &b, 00642 const std::vector<unsigned int> &variables, 00643 Number fptr(const Point& p, 00644 const Parameters& parameters, 00645 const std::string& sys_name, 00646 const std::string& unknown_name), 00647 Gradient gptr(const Point& p, 00648 const Parameters& parameters, 00649 const std::string& sys_name, 00650 const std::string& unknown_name), 00651 const Parameters& parameters) 00652 { 00653 WrappedFunction<Number> f(*this, fptr, ¶meters); 00654 WrappedFunction<Gradient> g(*this, gptr, ¶meters); 00655 this->boundary_project_solution(b, variables, &f, &g); 00656 } 00657 00658 00664 void System::boundary_project_solution 00665 (const std::set<boundary_id_type> &b, 00666 const std::vector<unsigned int> &variables, 00667 FunctionBase<Number> *f, 00668 FunctionBase<Gradient> *g) 00669 { 00670 this->boundary_project_vector(b, variables, *solution, f, g); 00671 00672 solution->localize(*current_local_solution); 00673 } 00674 00675 00676 00677 00678 00683 void System::boundary_project_vector 00684 (const std::set<boundary_id_type> &b, 00685 const std::vector<unsigned int> &variables, 00686 Number fptr(const Point& p, 00687 const Parameters& parameters, 00688 const std::string& sys_name, 00689 const std::string& unknown_name), 00690 Gradient gptr(const Point& p, 00691 const Parameters& parameters, 00692 const std::string& sys_name, 00693 const std::string& unknown_name), 00694 const Parameters& parameters, 00695 NumericVector<Number>& new_vector) const 00696 { 00697 WrappedFunction<Number> f(*this, fptr, ¶meters); 00698 WrappedFunction<Gradient> g(*this, gptr, ¶meters); 00699 this->boundary_project_vector(b, variables, new_vector, &f, &g); 00700 } 00701 00706 void System::boundary_project_vector 00707 (const std::set<boundary_id_type> &b, 00708 const std::vector<unsigned int> &variables, 00709 NumericVector<Number>& new_vector, 00710 FunctionBase<Number> *f, 00711 FunctionBase<Gradient> *g) const 00712 { 00713 START_LOG ("boundary_project_vector()", "System"); 00714 00715 Threads::parallel_for 00716 (ConstElemRange (this->get_mesh().active_local_elements_begin(), 00717 this->get_mesh().active_local_elements_end() ), 00718 BoundaryProjectSolution(b, variables, *this, f, g, 00719 this->get_equation_systems().parameters, 00720 new_vector) 00721 ); 00722 00723 // We don't do SCALAR dofs when just projecting the boundary, so 00724 // we're done here. 00725 00726 new_vector.close(); 00727 00728 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00729 this->get_dof_map().enforce_constraints_exactly(*this, &new_vector); 00730 #endif 00731 00732 STOP_LOG("boundary_project_vector()", "System"); 00733 } 00734 00735 00736 00737 #ifndef LIBMESH_ENABLE_AMR 00738 void ProjectVector::operator()(const ConstElemRange &) const 00739 { 00740 libmesh_error(); 00741 } 00742 #else 00743 void ProjectVector::operator()(const ConstElemRange &range) const 00744 { 00745 START_LOG ("operator()","ProjectVector"); 00746 00747 // A vector for Lagrange element interpolation, indicating if we 00748 // have visited a DOF yet. Note that this is thread-local storage, 00749 // hence shared DOFS that live on thread boundaries may be doubly 00750 // computed. It is expected that this will be more efficient 00751 // than locking a thread-global version of already_done, though. 00752 // 00753 // FIXME: we should use this for non-Lagrange coarsening, too 00754 std::vector<bool> already_done (system.n_dofs(), false); 00755 00756 // The number of variables in this system 00757 const unsigned int n_variables = system.n_vars(); 00758 00759 // The dimensionality of the current mesh 00760 const unsigned int dim = system.get_mesh().mesh_dimension(); 00761 00762 // The DofMap for this system 00763 const DofMap& dof_map = system.get_dof_map(); 00764 00765 // The element matrix and RHS for projections. 00766 // Note that Ke is always real-valued, whereas 00767 // Fe may be complex valued if complex number 00768 // support is enabled 00769 DenseMatrix<Real> Ke; 00770 DenseVector<Number> Fe; 00771 // The new element coefficients 00772 DenseVector<Number> Ue; 00773 00774 00775 // Loop over all the variables in the system 00776 for (unsigned int var=0; var<n_variables; var++) 00777 { 00778 const Variable& variable = dof_map.variable(var); 00779 00780 const FEType& base_fe_type = variable.type(); 00781 00782 if (base_fe_type.family == SCALAR) 00783 continue; 00784 00785 // Get FE objects of the appropriate type 00786 AutoPtr<FEBase> fe (FEBase::build(dim, base_fe_type)); 00787 AutoPtr<FEBase> fe_coarse (FEBase::build(dim, base_fe_type)); 00788 00789 // Create FE objects with potentially different p_level 00790 FEType fe_type, temp_fe_type; 00791 00792 // Prepare variables for non-Lagrange projection 00793 AutoPtr<QBase> qrule (base_fe_type.default_quadrature_rule(dim)); 00794 AutoPtr<QBase> qedgerule (base_fe_type.default_quadrature_rule(1)); 00795 AutoPtr<QBase> qsiderule (base_fe_type.default_quadrature_rule(dim-1)); 00796 std::vector<Point> coarse_qpoints; 00797 00798 // The values of the shape functions at the quadrature 00799 // points 00800 const std::vector<std::vector<Real> >& phi_values = 00801 fe->get_phi(); 00802 const std::vector<std::vector<Real> >& phi_coarse = 00803 fe_coarse->get_phi(); 00804 00805 // The Jacobian * quadrature weight at the quadrature points 00806 const std::vector<Real>& JxW = 00807 fe->get_JxW(); 00808 00809 // The XYZ locations of the quadrature points on the 00810 // child element 00811 const std::vector<Point>& xyz_values = 00812 fe->get_xyz(); 00813 00814 00815 // The global DOF indices 00816 std::vector<dof_id_type> new_dof_indices, old_dof_indices; 00817 // Side/edge local DOF indices 00818 std::vector<unsigned int> new_side_dofs, old_side_dofs; 00819 00820 // Iterate over the elements in the range 00821 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it) 00822 { 00823 const Elem* elem = *elem_it; 00824 // If this element doesn't have an old_dof_object with dofs for the 00825 // current system, then it must be newly added, so the user 00826 // is responsible for setting the new dofs. 00827 00828 // ... but we need a better way to test for that; the code 00829 // below breaks on any FE type for which the elem stores no 00830 // dofs. 00831 // if (!elem->old_dof_object || !elem->old_dof_object->has_dofs(system.number())) 00832 // continue; 00833 const Elem* parent = elem->parent(); 00834 00835 // Per-subdomain variables don't need to be projected on 00836 // elements where they're not active 00837 if (!variable.active_on_subdomain(elem->subdomain_id())) 00838 continue; 00839 00840 // Adjust the FE type for p-refined elements 00841 fe_type = base_fe_type; 00842 fe_type.order = static_cast<Order>(fe_type.order + 00843 elem->p_level()); 00844 00845 // We may need to remember the parent's p_level 00846 unsigned int old_parent_level = 0; 00847 00848 // Update the DOF indices for this element based on 00849 // the new mesh 00850 dof_map.dof_indices (elem, new_dof_indices, var); 00851 00852 // The number of DOFs on the new element 00853 const unsigned int new_n_dofs = 00854 libmesh_cast_int<unsigned int>(new_dof_indices.size()); 00855 00856 // Fixed vs. free DoFs on edge/face projections 00857 std::vector<char> dof_is_fixed(new_n_dofs, false); // bools 00858 std::vector<int> free_dof(new_n_dofs, 0); 00859 00860 // The element type 00861 const ElemType elem_type = elem->type(); 00862 00863 // The number of nodes on the new element 00864 const unsigned int n_nodes = elem->n_nodes(); 00865 00866 // Zero the interpolated values 00867 Ue.resize (new_n_dofs); Ue.zero(); 00868 00869 // Update the DOF indices based on the old mesh. 00870 // This is done in one of three ways: 00871 // 1.) If the child was just refined then it was not 00872 // active in the previous mesh & hence has no solution 00873 // values on it. In this case simply project or 00874 // interpolate the solution from the parent, who was 00875 // active in the previous mesh 00876 // 2.) If the child was just coarsened, obtaining a 00877 // well-defined solution may require doing independent 00878 // projections on nodes, edges, faces, and interiors 00879 // 3.) If the child was active in the previous 00880 // mesh, we can just copy coefficients directly 00881 if (elem->refinement_flag() == Elem::JUST_REFINED) 00882 { 00883 libmesh_assert(parent); 00884 old_parent_level = parent->p_level(); 00885 00886 // We may have done p refinement or coarsening as well; 00887 // if so then we need to reset the parent's p level 00888 // so we can get the right DoFs from it 00889 if (elem->p_refinement_flag() == Elem::JUST_REFINED) 00890 { 00891 libmesh_assert_greater (elem->p_level(), 0); 00892 (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() - 1); 00893 } 00894 else if (elem->p_refinement_flag() == Elem::JUST_COARSENED) 00895 { 00896 (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() + 1); 00897 } 00898 00899 dof_map.old_dof_indices (parent, old_dof_indices, var); 00900 } 00901 else 00902 { 00903 dof_map.old_dof_indices (elem, old_dof_indices, var); 00904 00905 if (elem->p_refinement_flag() == Elem::DO_NOTHING) 00906 libmesh_assert_equal_to (old_dof_indices.size(), new_n_dofs); 00907 if (elem->p_refinement_flag() == Elem::JUST_COARSENED) 00908 libmesh_assert (elem->has_children()); 00909 } 00910 00911 unsigned int old_n_dofs = 00912 libmesh_cast_int<unsigned int>(old_dof_indices.size()); 00913 00914 if (fe_type.family != LAGRANGE) { 00915 00916 // For refined non-Lagrange elements, we do an L2 00917 // projection 00918 // FIXME: this will be a suboptimal and ill-defined 00919 // result if we're using non-nested finite element 00920 // spaces or if we're on a p-coarsened element! 00921 if (elem->refinement_flag() == Elem::JUST_REFINED) 00922 { 00923 // Update the fe object based on the current child 00924 fe->attach_quadrature_rule (qrule.get()); 00925 fe->reinit (elem); 00926 00927 // The number of quadrature points on the child 00928 const unsigned int n_qp = qrule->n_points(); 00929 00930 FEInterface::inverse_map (dim, fe_type, parent, 00931 xyz_values, coarse_qpoints); 00932 00933 fe_coarse->reinit(parent, &coarse_qpoints); 00934 00935 // Reinitialize the element matrix and vector for 00936 // the current element. Note that this will zero them 00937 // before they are summed. 00938 Ke.resize (new_n_dofs, new_n_dofs); Ke.zero(); 00939 Fe.resize (new_n_dofs); Fe.zero(); 00940 00941 00942 // Loop over the quadrature points 00943 for (unsigned int qp=0; qp<n_qp; qp++) 00944 { 00945 // The solution value at the quadrature point 00946 Number val = libMesh::zero; 00947 00948 // Sum the function values * the DOF values 00949 // at the point of interest to get the function value 00950 // (Note that the # of DOFs on the parent need not be the 00951 // same as on the child!) 00952 for (unsigned int i=0; i<old_n_dofs; i++) 00953 { 00954 val += (old_vector(old_dof_indices[i])* 00955 phi_coarse[i][qp]); 00956 } 00957 00958 // Now \p val contains the solution value of variable 00959 // \p var at the qp'th quadrature point on the child 00960 // element. It has been interpolated from the parent 00961 // in case the child was just refined. Next we will 00962 // construct the L2-projection matrix for the element. 00963 00964 // Construct the Mass Matrix 00965 for (unsigned int i=0; i<new_n_dofs; i++) 00966 for (unsigned int j=0; j<new_n_dofs; j++) 00967 Ke(i,j) += JxW[qp]*phi_values[i][qp]*phi_values[j][qp]; 00968 00969 // Construct the RHS 00970 for (unsigned int i=0; i<new_n_dofs; i++) 00971 Fe(i) += JxW[qp]*phi_values[i][qp]*val; 00972 00973 } // end qp loop 00974 00975 Ke.cholesky_solve(Fe, Ue); 00976 00977 // Fix up the parent's p level in case we changed it 00978 (const_cast<Elem *>(parent))->hack_p_level(old_parent_level); 00979 } 00980 else if (elem->refinement_flag() == Elem::JUST_COARSENED) 00981 { 00982 FEBase::coarsened_dof_values(old_vector, dof_map, 00983 elem, Ue, var, true); 00984 } 00985 // For unrefined uncoarsened elements, we just copy DoFs 00986 else 00987 { 00988 // FIXME - I'm sure this function would be about half 00989 // the size if anyone ever figures out how to improve 00990 // the DofMap interface... - RHS 00991 if (elem->p_refinement_flag() == Elem::JUST_REFINED) 00992 { 00993 libmesh_assert_greater (elem->p_level(), 0); 00994 // P refinement of non-hierarchic bases will 00995 // require a whole separate code path 00996 libmesh_assert (fe->is_hierarchic()); 00997 temp_fe_type = fe_type; 00998 temp_fe_type.order = 00999 static_cast<Order>(temp_fe_type.order - 1); 01000 unsigned int old_index = 0, new_index = 0; 01001 for (unsigned int n=0; n != n_nodes; ++n) 01002 { 01003 const unsigned int nc = 01004 FEInterface::n_dofs_at_node (dim, temp_fe_type, 01005 elem_type, n); 01006 for (unsigned int i=0; i != nc; ++i) 01007 { 01008 Ue(new_index + i) = 01009 old_vector(old_dof_indices[old_index++]); 01010 } 01011 new_index += 01012 FEInterface::n_dofs_at_node (dim, fe_type, 01013 elem_type, n); 01014 } 01015 const unsigned int nc = 01016 FEInterface::n_dofs_per_elem (dim, temp_fe_type, 01017 elem_type); 01018 for (unsigned int i=0; i != nc; ++i) 01019 { 01020 Ue(new_index++) = 01021 old_vector(old_dof_indices[old_index+i]); 01022 } 01023 } 01024 else if (elem->p_refinement_flag() == 01025 Elem::JUST_COARSENED) 01026 { 01027 // P coarsening of non-hierarchic bases will 01028 // require a whole separate code path 01029 libmesh_assert (fe->is_hierarchic()); 01030 temp_fe_type = fe_type; 01031 temp_fe_type.order = 01032 static_cast<Order>(temp_fe_type.order + 1); 01033 unsigned int old_index = 0, new_index = 0; 01034 for (unsigned int n=0; n != n_nodes; ++n) 01035 { 01036 const unsigned int nc = 01037 FEInterface::n_dofs_at_node (dim, fe_type, 01038 elem_type, n); 01039 for (unsigned int i=0; i != nc; ++i) 01040 { 01041 Ue(new_index++) = 01042 old_vector(old_dof_indices[old_index+i]); 01043 } 01044 old_index += 01045 FEInterface::n_dofs_at_node (dim, temp_fe_type, 01046 elem_type, n); 01047 } 01048 const unsigned int nc = 01049 FEInterface::n_dofs_per_elem (dim, fe_type, 01050 elem_type); 01051 for (unsigned int i=0; i != nc; ++i) 01052 { 01053 Ue(new_index++) = 01054 old_vector(old_dof_indices[old_index+i]); 01055 } 01056 } 01057 else 01058 // If there's no p refinement, we can copy every DoF 01059 for (unsigned int i=0; i<new_n_dofs; i++) 01060 Ue(i) = old_vector(old_dof_indices[i]); 01061 } 01062 } 01063 else { // fe type is Lagrange 01064 // Loop over the DOFs on the element 01065 for (unsigned int new_local_dof=0; 01066 new_local_dof<new_n_dofs; new_local_dof++) 01067 { 01068 // The global DOF number for this local DOF 01069 const dof_id_type new_global_dof = 01070 new_dof_indices[new_local_dof]; 01071 01072 // The global DOF might lie outside of the bounds of a 01073 // distributed vector. Check for that and possibly continue 01074 // the loop 01075 if ((new_global_dof < new_vector.first_local_index()) || 01076 (new_global_dof >= new_vector.last_local_index())) 01077 continue; 01078 01079 // We might have already computed the solution for this DOF. 01080 // This is likely in the case of a shared node, particularly 01081 // at the corners of an element. Check to see if that is the 01082 // case 01083 if (already_done[new_global_dof] == true) 01084 continue; 01085 01086 already_done[new_global_dof] = true; 01087 01088 if (elem->refinement_flag() == Elem::JUST_REFINED) 01089 { 01090 // The location of the child's node on the parent element 01091 const Point point = 01092 FEInterface::inverse_map (dim, fe_type, parent, 01093 elem->point(new_local_dof)); 01094 01095 // Sum the function values * the DOF values 01096 // at the point of interest to get the function value 01097 // (Note that the # of DOFs on the parent need not be the 01098 // same as on the child!) 01099 for (unsigned int old_local_dof=0; 01100 old_local_dof<old_n_dofs; old_local_dof++) 01101 { 01102 const dof_id_type old_global_dof = 01103 old_dof_indices[old_local_dof]; 01104 01105 Ue(new_local_dof) += 01106 (old_vector(old_global_dof)* 01107 FEInterface::shape(dim, fe_type, parent, 01108 old_local_dof, point)); 01109 } 01110 } 01111 else 01112 { 01113 // Get the old global DOF index 01114 const dof_id_type old_global_dof = 01115 old_dof_indices[new_local_dof]; 01116 01117 Ue(new_local_dof) = old_vector(old_global_dof); 01118 } 01119 } // end local DOF loop 01120 01121 // We may have to clean up a parent's p_level 01122 if (elem->refinement_flag() == Elem::JUST_REFINED) 01123 (const_cast<Elem *>(parent))->hack_p_level(old_parent_level); 01124 } // end fe_type if() 01125 01126 // Lock the new_vector since it is shared among threads. 01127 { 01128 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 01129 01130 for (unsigned int i = 0; i < new_n_dofs; i++) 01131 if (Ue(i) != 0.) 01132 new_vector.set(new_dof_indices[i], Ue(i)); 01133 } 01134 } // end elem loop 01135 } // end variables loop 01136 01137 STOP_LOG ("operator()","ProjectVector"); 01138 } 01139 #endif // LIBMESH_ENABLE_AMR 01140 01141 01142 01143 void BuildProjectionList::unique() 01144 { 01145 // Sort the send list. After this duplicated 01146 // elements will be adjacent in the vector 01147 std::sort(this->send_list.begin(), 01148 this->send_list.end()); 01149 01150 // Now use std::unique to remove duplicate entries 01151 std::vector<dof_id_type>::iterator new_end = 01152 std::unique (this->send_list.begin(), 01153 this->send_list.end()); 01154 01155 // Remove the end of the send_list. Use the "swap trick" 01156 // from Effective STL 01157 std::vector<dof_id_type> 01158 (this->send_list.begin(), new_end).swap (this->send_list); 01159 } 01160 01161 01162 01163 #ifndef LIBMESH_ENABLE_AMR 01164 void BuildProjectionList::operator()(const ConstElemRange &) 01165 { 01166 libmesh_error(); 01167 } 01168 #else 01169 void BuildProjectionList::operator()(const ConstElemRange &range) 01170 { 01171 // The DofMap for this system 01172 const DofMap& dof_map = system.get_dof_map(); 01173 01174 const dof_id_type first_old_dof = dof_map.first_old_dof(); 01175 const dof_id_type end_old_dof = dof_map.end_old_dof(); 01176 01177 // We can handle all the variables at once. 01178 // The old global DOF indices 01179 std::vector<dof_id_type> di; 01180 01181 // Iterate over the elements in the range 01182 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it) 01183 { 01184 const Elem* elem = *elem_it; 01185 // If this element doesn't have an old_dof_object with dofs for the 01186 // current system, then it must be newly added, so the user 01187 // is responsible for setting the new dofs. 01188 01189 // ... but we need a better way to test for that; the code 01190 // below breaks on any FE type for which the elem stores no 01191 // dofs. 01192 // if (!elem->old_dof_object || !elem->old_dof_object->has_dofs(system.number())) 01193 // continue; 01194 const Elem* parent = elem->parent(); 01195 01196 if (elem->refinement_flag() == Elem::JUST_REFINED) 01197 { 01198 libmesh_assert(parent); 01199 unsigned int old_parent_level = parent->p_level(); 01200 01201 if (elem->p_refinement_flag() == Elem::JUST_REFINED) 01202 { 01203 // We may have done p refinement or coarsening as well; 01204 // if so then we need to reset the parent's p level 01205 // so we can get the right DoFs from it 01206 libmesh_assert_greater (elem->p_level(), 0); 01207 (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() - 1); 01208 } 01209 else if (elem->p_refinement_flag() == Elem::JUST_COARSENED) 01210 { 01211 (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() + 1); 01212 } 01213 01214 dof_map.old_dof_indices (parent, di); 01215 01216 // Fix up the parent's p level in case we changed it 01217 (const_cast<Elem *>(parent))->hack_p_level(old_parent_level); 01218 } 01219 else if (elem->refinement_flag() == Elem::JUST_COARSENED) 01220 { 01221 std::vector<dof_id_type> di_child; 01222 di.clear(); 01223 for (unsigned int c=0; c != elem->n_children(); ++c) 01224 { 01225 dof_map.old_dof_indices (elem->child(c), di_child); 01226 di.insert(di.end(), di_child.begin(), di_child.end()); 01227 } 01228 } 01229 else 01230 dof_map.old_dof_indices (elem, di); 01231 01232 for (unsigned int i=0; i != di.size(); ++i) 01233 if (di[i] < first_old_dof || di[i] >= end_old_dof) 01234 this->send_list.push_back(di[i]); 01235 } // end elem loop 01236 } 01237 #endif // LIBMESH_ENABLE_AMR 01238 01239 01240 01241 void BuildProjectionList::join(const BuildProjectionList &other) 01242 { 01243 // Joining simply requires I add the dof indices from the other object 01244 this->send_list.insert(this->send_list.end(), 01245 other.send_list.begin(), 01246 other.send_list.end()); 01247 } 01248 01249 01250 void ProjectSolution::operator()(const ConstElemRange &range) const 01251 { 01252 // We need data to project 01253 libmesh_assert(f.get()); 01254 01262 // The number of variables in this system 01263 const unsigned int n_variables = system.n_vars(); 01264 01265 // The dimensionality of the current mesh 01266 const unsigned int dim = system.get_mesh().mesh_dimension(); 01267 01268 // The DofMap for this system 01269 const DofMap& dof_map = system.get_dof_map(); 01270 01271 // The element matrix and RHS for projections. 01272 // Note that Ke is always real-valued, whereas 01273 // Fe may be complex valued if complex number 01274 // support is enabled 01275 DenseMatrix<Real> Ke; 01276 DenseVector<Number> Fe; 01277 // The new element coefficients 01278 DenseVector<Number> Ue; 01279 01280 01281 // Loop over all the variables in the system 01282 for (unsigned int var=0; var<n_variables; var++) 01283 { 01284 const Variable& variable = dof_map.variable(var); 01285 01286 const FEType& fe_type = variable.type(); 01287 01288 if (fe_type.family == SCALAR) 01289 continue; 01290 01291 const unsigned int var_component = 01292 system.variable_scalar_number(var, 0); 01293 01294 // Get FE objects of the appropriate type 01295 AutoPtr<FEBase> fe (FEBase::build(dim, fe_type)); 01296 01297 // Prepare variables for projection 01298 AutoPtr<QBase> qrule (fe_type.default_quadrature_rule(dim)); 01299 AutoPtr<QBase> qedgerule (fe_type.default_quadrature_rule(1)); 01300 AutoPtr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1)); 01301 01302 // The values of the shape functions at the quadrature 01303 // points 01304 const std::vector<std::vector<Real> >& phi = fe->get_phi(); 01305 01306 // The gradients of the shape functions at the quadrature 01307 // points on the child element. 01308 const std::vector<std::vector<RealGradient> > *dphi = NULL; 01309 01310 const FEContinuity cont = fe->get_continuity(); 01311 01312 if (cont == C_ONE) 01313 { 01314 // We'll need gradient data for a C1 projection 01315 libmesh_assert(g.get()); 01316 01317 const std::vector<std::vector<RealGradient> >& 01318 ref_dphi = fe->get_dphi(); 01319 dphi = &ref_dphi; 01320 } 01321 01322 // The Jacobian * quadrature weight at the quadrature points 01323 const std::vector<Real>& JxW = 01324 fe->get_JxW(); 01325 01326 // The XYZ locations of the quadrature points 01327 const std::vector<Point>& xyz_values = 01328 fe->get_xyz(); 01329 01330 // The global DOF indices 01331 std::vector<dof_id_type> dof_indices; 01332 // Side/edge DOF indices 01333 std::vector<unsigned int> side_dofs; 01334 01335 // Iterate over all the elements in the range 01336 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it) 01337 { 01338 const Elem* elem = *elem_it; 01339 01340 // Per-subdomain variables don't need to be projected on 01341 // elements where they're not active 01342 if (!variable.active_on_subdomain(elem->subdomain_id())) 01343 continue; 01344 01345 // Update the DOF indices for this element based on 01346 // the current mesh 01347 dof_map.dof_indices (elem, dof_indices, var); 01348 01349 // The number of DOFs on the element 01350 const unsigned int n_dofs = 01351 libmesh_cast_int<unsigned int>(dof_indices.size()); 01352 01353 // Fixed vs. free DoFs on edge/face projections 01354 std::vector<char> dof_is_fixed(n_dofs, false); // bools 01355 std::vector<int> free_dof(n_dofs, 0); 01356 01357 // The element type 01358 const ElemType elem_type = elem->type(); 01359 01360 // The number of nodes on the new element 01361 const unsigned int n_nodes = elem->n_nodes(); 01362 01363 // Zero the interpolated values 01364 Ue.resize (n_dofs); Ue.zero(); 01365 01366 // In general, we need a series of 01367 // projections to ensure a unique and continuous 01368 // solution. We start by interpolating nodes, then 01369 // hold those fixed and project edges, then 01370 // hold those fixed and project faces, then 01371 // hold those fixed and project interiors 01372 01373 // Interpolate node values first 01374 unsigned int current_dof = 0; 01375 for (unsigned int n=0; n!= n_nodes; ++n) 01376 { 01377 // FIXME: this should go through the DofMap, 01378 // not duplicate dof_indices code badly! 01379 const unsigned int nc = 01380 FEInterface::n_dofs_at_node (dim, fe_type, elem_type, 01381 n); 01382 if (!elem->is_vertex(n)) 01383 { 01384 current_dof += nc; 01385 continue; 01386 } 01387 if (cont == DISCONTINUOUS) 01388 { 01389 libmesh_assert_equal_to (nc, 0); 01390 } 01391 // Assume that C_ZERO elements have a single nodal 01392 // value shape function 01393 else if (cont == C_ZERO) 01394 { 01395 libmesh_assert_equal_to (nc, 1); 01396 Ue(current_dof) = f->component(var_component, 01397 elem->point(n), 01398 system.time); 01399 dof_is_fixed[current_dof] = true; 01400 current_dof++; 01401 } 01402 // The hermite element vertex shape functions are weird 01403 else if (fe_type.family == HERMITE) 01404 { 01405 Ue(current_dof) = f->component(var_component, 01406 elem->point(n), 01407 system.time); 01408 dof_is_fixed[current_dof] = true; 01409 current_dof++; 01410 Gradient grad = g->component(var_component, 01411 elem->point(n), 01412 system.time); 01413 // x derivative 01414 Ue(current_dof) = grad(0); 01415 dof_is_fixed[current_dof] = true; 01416 current_dof++; 01417 if (dim > 1) 01418 { 01419 // We'll finite difference mixed derivatives 01420 Point nxminus = elem->point(n), 01421 nxplus = elem->point(n); 01422 nxminus(0) -= TOLERANCE; 01423 nxplus(0) += TOLERANCE; 01424 Gradient gxminus = g->component(var_component, 01425 nxminus, 01426 system.time); 01427 Gradient gxplus = g->component(var_component, 01428 nxplus, 01429 system.time); 01430 // y derivative 01431 Ue(current_dof) = grad(1); 01432 dof_is_fixed[current_dof] = true; 01433 current_dof++; 01434 // xy derivative 01435 Ue(current_dof) = (gxplus(1) - gxminus(1)) 01436 / 2. / TOLERANCE; 01437 dof_is_fixed[current_dof] = true; 01438 current_dof++; 01439 01440 if (dim > 2) 01441 { 01442 // z derivative 01443 Ue(current_dof) = grad(2); 01444 dof_is_fixed[current_dof] = true; 01445 current_dof++; 01446 // xz derivative 01447 Ue(current_dof) = (gxplus(2) - gxminus(2)) 01448 / 2. / TOLERANCE; 01449 dof_is_fixed[current_dof] = true; 01450 current_dof++; 01451 // We need new points for yz 01452 Point nyminus = elem->point(n), 01453 nyplus = elem->point(n); 01454 nyminus(1) -= TOLERANCE; 01455 nyplus(1) += TOLERANCE; 01456 Gradient gyminus = g->component(var_component, 01457 nyminus, 01458 system.time); 01459 Gradient gyplus = g->component(var_component, 01460 nyplus, 01461 system.time); 01462 // xz derivative 01463 Ue(current_dof) = (gyplus(2) - gyminus(2)) 01464 / 2. / TOLERANCE; 01465 dof_is_fixed[current_dof] = true; 01466 current_dof++; 01467 // Getting a 2nd order xyz is more tedious 01468 Point nxmym = elem->point(n), 01469 nxmyp = elem->point(n), 01470 nxpym = elem->point(n), 01471 nxpyp = elem->point(n); 01472 nxmym(0) -= TOLERANCE; 01473 nxmym(1) -= TOLERANCE; 01474 nxmyp(0) -= TOLERANCE; 01475 nxmyp(1) += TOLERANCE; 01476 nxpym(0) += TOLERANCE; 01477 nxpym(1) -= TOLERANCE; 01478 nxpyp(0) += TOLERANCE; 01479 nxpyp(1) += TOLERANCE; 01480 Gradient gxmym = g->component(var_component, 01481 nxmym, 01482 system.time); 01483 Gradient gxmyp = g->component(var_component, 01484 nxmyp, 01485 system.time); 01486 Gradient gxpym = g->component(var_component, 01487 nxpym, 01488 system.time); 01489 Gradient gxpyp = g->component(var_component, 01490 nxpyp, 01491 system.time); 01492 Number gxzplus = (gxpyp(2) - gxmyp(2)) 01493 / 2. / TOLERANCE; 01494 Number gxzminus = (gxpym(2) - gxmym(2)) 01495 / 2. / TOLERANCE; 01496 // xyz derivative 01497 Ue(current_dof) = (gxzplus - gxzminus) 01498 / 2. / TOLERANCE; 01499 dof_is_fixed[current_dof] = true; 01500 current_dof++; 01501 } 01502 } 01503 } 01504 // Assume that other C_ONE elements have a single nodal 01505 // value shape function and nodal gradient component 01506 // shape functions 01507 else if (cont == C_ONE) 01508 { 01509 libmesh_assert_equal_to (nc, 1 + dim); 01510 Ue(current_dof) = f->component(var_component, 01511 elem->point(n), 01512 system.time); 01513 dof_is_fixed[current_dof] = true; 01514 current_dof++; 01515 Gradient grad = g->component(var_component, 01516 elem->point(n), 01517 system.time); 01518 for (unsigned int i=0; i!= dim; ++i) 01519 { 01520 Ue(current_dof) = grad(i); 01521 dof_is_fixed[current_dof] = true; 01522 current_dof++; 01523 } 01524 } 01525 else 01526 libmesh_error(); 01527 } 01528 01529 // In 3D, project any edge values next 01530 if (dim > 2 && cont != DISCONTINUOUS) 01531 for (unsigned int e=0; e != elem->n_edges(); ++e) 01532 { 01533 FEInterface::dofs_on_edge(elem, dim, fe_type, e, 01534 side_dofs); 01535 01536 // Some edge dofs are on nodes and already 01537 // fixed, others are free to calculate 01538 unsigned int free_dofs = 0; 01539 for (unsigned int i=0; i != side_dofs.size(); ++i) 01540 if (!dof_is_fixed[side_dofs[i]]) 01541 free_dof[free_dofs++] = i; 01542 01543 // There may be nothing to project 01544 if (!free_dofs) 01545 continue; 01546 01547 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01548 Fe.resize (free_dofs); Fe.zero(); 01549 // The new edge coefficients 01550 DenseVector<Number> Uedge(free_dofs); 01551 01552 // Initialize FE data on the edge 01553 fe->attach_quadrature_rule (qedgerule.get()); 01554 fe->edge_reinit (elem, e); 01555 const unsigned int n_qp = qedgerule->n_points(); 01556 01557 // Loop over the quadrature points 01558 for (unsigned int qp=0; qp<n_qp; qp++) 01559 { 01560 // solution at the quadrature point 01561 Number fineval = f->component(var_component, 01562 xyz_values[qp], 01563 system.time); 01564 // solution grad at the quadrature point 01565 Gradient finegrad; 01566 if (cont == C_ONE) 01567 finegrad = g->component(var_component, 01568 xyz_values[qp], 01569 system.time); 01570 01571 // Form edge projection matrix 01572 for (unsigned int sidei=0, freei=0; 01573 sidei != side_dofs.size(); ++sidei) 01574 { 01575 unsigned int i = side_dofs[sidei]; 01576 // fixed DoFs aren't test functions 01577 if (dof_is_fixed[i]) 01578 continue; 01579 for (unsigned int sidej=0, freej=0; 01580 sidej != side_dofs.size(); ++sidej) 01581 { 01582 unsigned int j = side_dofs[sidej]; 01583 if (dof_is_fixed[j]) 01584 Fe(freei) -= phi[i][qp] * phi[j][qp] * 01585 JxW[qp] * Ue(j); 01586 else 01587 Ke(freei,freej) += phi[i][qp] * 01588 phi[j][qp] * JxW[qp]; 01589 if (cont == C_ONE) 01590 { 01591 if (dof_is_fixed[j]) 01592 Fe(freei) -= ((*dphi)[i][qp] * 01593 (*dphi)[j][qp]) * 01594 JxW[qp] * Ue(j); 01595 else 01596 Ke(freei,freej) += ((*dphi)[i][qp] * 01597 (*dphi)[j][qp]) 01598 * JxW[qp]; 01599 } 01600 if (!dof_is_fixed[j]) 01601 freej++; 01602 } 01603 Fe(freei) += phi[i][qp] * fineval * JxW[qp]; 01604 if (cont == C_ONE) 01605 Fe(freei) += (finegrad * (*dphi)[i][qp]) * 01606 JxW[qp]; 01607 freei++; 01608 } 01609 } 01610 01611 Ke.cholesky_solve(Fe, Uedge); 01612 01613 // Transfer new edge solutions to element 01614 for (unsigned int i=0; i != free_dofs; ++i) 01615 { 01616 Number &ui = Ue(side_dofs[free_dof[i]]); 01617 libmesh_assert(std::abs(ui) < TOLERANCE || 01618 std::abs(ui - Uedge(i)) < TOLERANCE); 01619 ui = Uedge(i); 01620 dof_is_fixed[side_dofs[free_dof[i]]] = true; 01621 } 01622 } 01623 01624 // Project any side values (edges in 2D, faces in 3D) 01625 if (dim > 1 && cont != DISCONTINUOUS) 01626 for (unsigned int s=0; s != elem->n_sides(); ++s) 01627 { 01628 FEInterface::dofs_on_side(elem, dim, fe_type, s, 01629 side_dofs); 01630 01631 // Some side dofs are on nodes/edges and already 01632 // fixed, others are free to calculate 01633 unsigned int free_dofs = 0; 01634 for (unsigned int i=0; i != side_dofs.size(); ++i) 01635 if (!dof_is_fixed[side_dofs[i]]) 01636 free_dof[free_dofs++] = i; 01637 01638 // There may be nothing to project 01639 if (!free_dofs) 01640 continue; 01641 01642 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01643 Fe.resize (free_dofs); Fe.zero(); 01644 // The new side coefficients 01645 DenseVector<Number> Uside(free_dofs); 01646 01647 // Initialize FE data on the side 01648 fe->attach_quadrature_rule (qsiderule.get()); 01649 fe->reinit (elem, s); 01650 const unsigned int n_qp = qsiderule->n_points(); 01651 01652 // Loop over the quadrature points 01653 for (unsigned int qp=0; qp<n_qp; qp++) 01654 { 01655 // solution at the quadrature point 01656 Number fineval = f->component(var_component, 01657 xyz_values[qp], 01658 system.time); 01659 // solution grad at the quadrature point 01660 Gradient finegrad; 01661 if (cont == C_ONE) 01662 finegrad = g->component(var_component, 01663 xyz_values[qp], 01664 system.time); 01665 01666 // Form side projection matrix 01667 for (unsigned int sidei=0, freei=0; 01668 sidei != side_dofs.size(); ++sidei) 01669 { 01670 unsigned int i = side_dofs[sidei]; 01671 // fixed DoFs aren't test functions 01672 if (dof_is_fixed[i]) 01673 continue; 01674 for (unsigned int sidej=0, freej=0; 01675 sidej != side_dofs.size(); ++sidej) 01676 { 01677 unsigned int j = side_dofs[sidej]; 01678 if (dof_is_fixed[j]) 01679 Fe(freei) -= phi[i][qp] * phi[j][qp] * 01680 JxW[qp] * Ue(j); 01681 else 01682 Ke(freei,freej) += phi[i][qp] * 01683 phi[j][qp] * JxW[qp]; 01684 if (cont == C_ONE) 01685 { 01686 if (dof_is_fixed[j]) 01687 Fe(freei) -= ((*dphi)[i][qp] * 01688 (*dphi)[j][qp]) * 01689 JxW[qp] * Ue(j); 01690 else 01691 Ke(freei,freej) += ((*dphi)[i][qp] * 01692 (*dphi)[j][qp]) 01693 * JxW[qp]; 01694 } 01695 if (!dof_is_fixed[j]) 01696 freej++; 01697 } 01698 Fe(freei) += (fineval * phi[i][qp]) * JxW[qp]; 01699 if (cont == C_ONE) 01700 Fe(freei) += (finegrad * (*dphi)[i][qp]) * 01701 JxW[qp]; 01702 freei++; 01703 } 01704 } 01705 01706 Ke.cholesky_solve(Fe, Uside); 01707 01708 // Transfer new side solutions to element 01709 for (unsigned int i=0; i != free_dofs; ++i) 01710 { 01711 Number &ui = Ue(side_dofs[free_dof[i]]); 01712 libmesh_assert(std::abs(ui) < TOLERANCE || 01713 std::abs(ui - Uside(i)) < TOLERANCE); 01714 ui = Uside(i); 01715 dof_is_fixed[side_dofs[free_dof[i]]] = true; 01716 } 01717 } 01718 01719 // Project the interior values, finally 01720 01721 // Some interior dofs are on nodes/edges/sides and 01722 // already fixed, others are free to calculate 01723 unsigned int free_dofs = 0; 01724 for (unsigned int i=0; i != n_dofs; ++i) 01725 if (!dof_is_fixed[i]) 01726 free_dof[free_dofs++] = i; 01727 01728 // There may be nothing to project 01729 if (free_dofs) 01730 { 01731 01732 Ke.resize (free_dofs, free_dofs); Ke.zero(); 01733 Fe.resize (free_dofs); Fe.zero(); 01734 // The new interior coefficients 01735 DenseVector<Number> Uint(free_dofs); 01736 01737 // Initialize FE data 01738 fe->attach_quadrature_rule (qrule.get()); 01739 fe->reinit (elem); 01740 const unsigned int n_qp = qrule->n_points(); 01741 01742 // Loop over the quadrature points 01743 for (unsigned int qp=0; qp<n_qp; qp++) 01744 { 01745 // solution at the quadrature point 01746 Number fineval = f->component(var_component, 01747 xyz_values[qp], 01748 system.time); 01749 // solution grad at the quadrature point 01750 Gradient finegrad; 01751 if (cont == C_ONE) 01752 finegrad = g->component(var_component, 01753 xyz_values[qp], 01754 system.time); 01755 01756 // Form interior projection matrix 01757 for (unsigned int i=0, freei=0; i != n_dofs; ++i) 01758 { 01759 // fixed DoFs aren't test functions 01760 if (dof_is_fixed[i]) 01761 continue; 01762 for (unsigned int j=0, freej=0; j != n_dofs; ++j) 01763 { 01764 if (dof_is_fixed[j]) 01765 Fe(freei) -= phi[i][qp] * phi[j][qp] * JxW[qp] 01766 * Ue(j); 01767 else 01768 Ke(freei,freej) += phi[i][qp] * phi[j][qp] * 01769 JxW[qp]; 01770 if (cont == C_ONE) 01771 { 01772 if (dof_is_fixed[j]) 01773 Fe(freei) -= ((*dphi)[i][qp] * 01774 (*dphi)[j][qp]) * JxW[qp] * 01775 Ue(j); 01776 else 01777 Ke(freei,freej) += ((*dphi)[i][qp] * 01778 (*dphi)[j][qp]) * 01779 JxW[qp]; 01780 } 01781 if (!dof_is_fixed[j]) 01782 freej++; 01783 } 01784 Fe(freei) += phi[i][qp] * fineval * JxW[qp]; 01785 if (cont == C_ONE) 01786 Fe(freei) += (finegrad * (*dphi)[i][qp]) * JxW[qp]; 01787 freei++; 01788 } 01789 } 01790 Ke.cholesky_solve(Fe, Uint); 01791 01792 // Transfer new interior solutions to element 01793 for (unsigned int i=0; i != free_dofs; ++i) 01794 { 01795 Number &ui = Ue(free_dof[i]); 01796 libmesh_assert(std::abs(ui) < TOLERANCE || 01797 std::abs(ui - Uint(i)) < TOLERANCE); 01798 ui = Uint(i); 01799 dof_is_fixed[free_dof[i]] = true; 01800 } 01801 01802 } // if there are free interior dofs 01803 01804 // Make sure every DoF got reached! 01805 for (unsigned int i=0; i != n_dofs; ++i) 01806 libmesh_assert(dof_is_fixed[i]); 01807 01808 const dof_id_type 01809 first = new_vector.first_local_index(), 01810 last = new_vector.last_local_index(); 01811 01812 // Lock the new_vector since it is shared among threads. 01813 { 01814 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 01815 01816 for (unsigned int i = 0; i < n_dofs; i++) 01817 // We may be projecting a new zero value onto 01818 // an old nonzero approximation - RHS 01819 // if (Ue(i) != 0.) 01820 if ((dof_indices[i] >= first) && 01821 (dof_indices[i] < last)) 01822 { 01823 new_vector.set(dof_indices[i], Ue(i)); 01824 } 01825 } 01826 } // end elem loop 01827 } // end variables loop 01828 } 01829 01830 01831 void ProjectFEMSolution::operator()(const ConstElemRange &range) const 01832 { 01833 // We need data to project 01834 libmesh_assert(f.get()); 01835 01843 FEMContext context( system ); 01844 01845 // The number of variables in this system 01846 const unsigned int n_variables = context.n_vars(); 01847 01848 // The dimensionality of the current mesh 01849 const unsigned int dim = context.get_dim(); 01850 01851 // The DofMap for this system 01852 const DofMap& dof_map = system.get_dof_map(); 01853 01854 // The element matrix and RHS for projections. 01855 // Note that Ke is always real-valued, whereas 01856 // Fe may be complex valued if complex number 01857 // support is enabled 01858 DenseMatrix<Real> Ke; 01859 DenseVector<Number> Fe; 01860 // The new element coefficients 01861 DenseVector<Number> Ue; 01862 01863 // FIXME: Need to generalize this to vector-valued elements. [PB] 01864 FEBase* fe = NULL; 01865 FEBase* side_fe = NULL; 01866 FEBase* edge_fe = NULL; 01867 01868 // First, loop over all variables and make sure the shape functions phi will 01869 // get built as well as the shape function gradients if the gradient functor 01870 // is supplied. 01871 for (unsigned int var=0; var<n_variables; var++) 01872 { 01873 context.get_element_fe( var, fe ); 01874 if (fe->get_fe_type().family == SCALAR) 01875 continue; 01876 if( dim > 1 ) 01877 context.get_side_fe( var, side_fe ); 01878 if( dim > 2 ) 01879 context.get_edge_fe( var, edge_fe ); 01880 01881 fe->get_phi(); 01882 if( dim > 1 ) 01883 side_fe->get_phi(); 01884 if( dim > 2 ) 01885 edge_fe->get_phi(); 01886 01887 const FEContinuity cont = fe->get_continuity(); 01888 if(cont == C_ONE) 01889 { 01890 libmesh_assert(g.get()); 01891 fe->get_dphi(); 01892 side_fe->get_dphi(); 01893 if( dim > 2 ) 01894 edge_fe->get_dphi(); 01895 } 01896 } 01897 01898 // Now initialize any user requested shape functions 01899 f->init_context(context); 01900 if(g.get()) 01901 g->init_context(context); 01902 01903 std::vector<unsigned int> side_dofs; 01904 01905 // Iterate over all the elements in the range 01906 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it) 01907 { 01908 const Elem* elem = *elem_it; 01909 01910 context.pre_fe_reinit(system, elem); 01911 01912 // Loop over all the variables in the system 01913 for (unsigned int var=0; var<n_variables; var++) 01914 { 01915 const Variable& variable = dof_map.variable(var); 01916 01917 const FEType& fe_type = variable.type(); 01918 01919 if (fe_type.family == SCALAR) 01920 continue; 01921 01922 // Per-subdomain variables don't need to be projected on 01923 // elements where they're not active 01924 if (!variable.active_on_subdomain(elem->subdomain_id())) 01925 continue; 01926 01927 const FEContinuity cont = fe->get_continuity(); 01928 01929 const unsigned int var_component = 01930 system.variable_scalar_number(var, 0); 01931 01932 const std::vector<dof_id_type>& dof_indices = 01933 context.get_dof_indices(var); 01934 01935 // The number of DOFs on the element 01936 const unsigned int n_dofs = 01937 libmesh_cast_int<unsigned int>(dof_indices.size()); 01938 01939 // Fixed vs. free DoFs on edge/face projections 01940 std::vector<char> dof_is_fixed(n_dofs, false); // bools 01941 std::vector<int> free_dof(n_dofs, 0); 01942 01943 // The element type 01944 const ElemType elem_type = elem->type(); 01945 01946 // The number of nodes on the new element 01947 const unsigned int n_nodes = elem->n_nodes(); 01948 01949 // Zero the interpolated values 01950 Ue.resize (n_dofs); Ue.zero(); 01951 01952 // In general, we need a series of 01953 // projections to ensure a unique and continuous 01954 // solution. We start by interpolating nodes, then 01955 // hold those fixed and project edges, then 01956 // hold those fixed and project faces, then 01957 // hold those fixed and project interiors 01958 01959 // Interpolate node values first 01960 unsigned int current_dof = 0; 01961 for (unsigned int n=0; n!= n_nodes; ++n) 01962 { 01963 // FIXME: this should go through the DofMap, 01964 // not duplicate dof_indices code badly! 01965 const unsigned int nc = 01966 FEInterface::n_dofs_at_node (dim, fe_type, elem_type, 01967 n); 01968 if (!elem->is_vertex(n)) 01969 { 01970 current_dof += nc; 01971 continue; 01972 } 01973 if (cont == DISCONTINUOUS) 01974 { 01975 libmesh_assert_equal_to (nc, 0); 01976 } 01977 // Assume that C_ZERO elements have a single nodal 01978 // value shape function 01979 else if (cont == C_ZERO) 01980 { 01981 libmesh_assert_equal_to (nc, 1); 01982 Ue(current_dof) = f->component(context, 01983 var_component, 01984 elem->point(n), 01985 system.time); 01986 dof_is_fixed[current_dof] = true; 01987 current_dof++; 01988 } 01989 // The hermite element vertex shape functions are weird 01990 else if (fe_type.family == HERMITE) 01991 { 01992 Ue(current_dof) = f->component(context, 01993 var_component, 01994 elem->point(n), 01995 system.time); 01996 dof_is_fixed[current_dof] = true; 01997 current_dof++; 01998 Gradient grad = g->component(context, 01999 var_component, 02000 elem->point(n), 02001 system.time); 02002 // x derivative 02003 Ue(current_dof) = grad(0); 02004 dof_is_fixed[current_dof] = true; 02005 current_dof++; 02006 if (dim > 1) 02007 { 02008 // We'll finite difference mixed derivatives 02009 Point nxminus = elem->point(n), 02010 nxplus = elem->point(n); 02011 nxminus(0) -= TOLERANCE; 02012 nxplus(0) += TOLERANCE; 02013 Gradient gxminus = g->component(context, 02014 var_component, 02015 nxminus, 02016 system.time); 02017 Gradient gxplus = g->component(context, 02018 var_component, 02019 nxplus, 02020 system.time); 02021 // y derivative 02022 Ue(current_dof) = grad(1); 02023 dof_is_fixed[current_dof] = true; 02024 current_dof++; 02025 // xy derivative 02026 Ue(current_dof) = (gxplus(1) - gxminus(1)) 02027 / 2. / TOLERANCE; 02028 dof_is_fixed[current_dof] = true; 02029 current_dof++; 02030 02031 if (dim > 2) 02032 { 02033 // z derivative 02034 Ue(current_dof) = grad(2); 02035 dof_is_fixed[current_dof] = true; 02036 current_dof++; 02037 // xz derivative 02038 Ue(current_dof) = (gxplus(2) - gxminus(2)) 02039 / 2. / TOLERANCE; 02040 dof_is_fixed[current_dof] = true; 02041 current_dof++; 02042 // We need new points for yz 02043 Point nyminus = elem->point(n), 02044 nyplus = elem->point(n); 02045 nyminus(1) -= TOLERANCE; 02046 nyplus(1) += TOLERANCE; 02047 Gradient gyminus = g->component(context, 02048 var_component, 02049 nyminus, 02050 system.time); 02051 Gradient gyplus = g->component(context, 02052 var_component, 02053 nyplus, 02054 system.time); 02055 // xz derivative 02056 Ue(current_dof) = (gyplus(2) - gyminus(2)) 02057 / 2. / TOLERANCE; 02058 dof_is_fixed[current_dof] = true; 02059 current_dof++; 02060 // Getting a 2nd order xyz is more tedious 02061 Point nxmym = elem->point(n), 02062 nxmyp = elem->point(n), 02063 nxpym = elem->point(n), 02064 nxpyp = elem->point(n); 02065 nxmym(0) -= TOLERANCE; 02066 nxmym(1) -= TOLERANCE; 02067 nxmyp(0) -= TOLERANCE; 02068 nxmyp(1) += TOLERANCE; 02069 nxpym(0) += TOLERANCE; 02070 nxpym(1) -= TOLERANCE; 02071 nxpyp(0) += TOLERANCE; 02072 nxpyp(1) += TOLERANCE; 02073 Gradient gxmym = g->component(context, 02074 var_component, 02075 nxmym, 02076 system.time); 02077 Gradient gxmyp = g->component(context, 02078 var_component, 02079 nxmyp, 02080 system.time); 02081 Gradient gxpym = g->component(context, 02082 var_component, 02083 nxpym, 02084 system.time); 02085 Gradient gxpyp = g->component(context, 02086 var_component, 02087 nxpyp, 02088 system.time); 02089 Number gxzplus = (gxpyp(2) - gxmyp(2)) 02090 / 2. / TOLERANCE; 02091 Number gxzminus = (gxpym(2) - gxmym(2)) 02092 / 2. / TOLERANCE; 02093 // xyz derivative 02094 Ue(current_dof) = (gxzplus - gxzminus) 02095 / 2. / TOLERANCE; 02096 dof_is_fixed[current_dof] = true; 02097 current_dof++; 02098 } 02099 } 02100 } 02101 // Assume that other C_ONE elements have a single nodal 02102 // value shape function and nodal gradient component 02103 // shape functions 02104 else if (cont == C_ONE) 02105 { 02106 libmesh_assert_equal_to (nc, 1 + dim); 02107 Ue(current_dof) = f->component(context, 02108 var_component, 02109 elem->point(n), 02110 system.time); 02111 dof_is_fixed[current_dof] = true; 02112 current_dof++; 02113 Gradient grad = g->component(context, 02114 var_component, 02115 elem->point(n), 02116 system.time); 02117 for (unsigned int i=0; i!= dim; ++i) 02118 { 02119 Ue(current_dof) = grad(i); 02120 dof_is_fixed[current_dof] = true; 02121 current_dof++; 02122 } 02123 } 02124 else 02125 libmesh_error(); 02126 } 02127 02128 // In 3D, project any edge values next 02129 if (dim > 2 && cont != DISCONTINUOUS) 02130 { 02131 const QBase* qedgerule = context.get_edge_qrule(); 02132 const unsigned int n_qp = qedgerule->n_points(); 02133 const std::vector<Point>& xyz_values = edge_fe->get_xyz(); 02134 const std::vector<Real>& JxW = edge_fe->get_JxW(); 02135 02136 const std::vector<std::vector<Real> >& phi = edge_fe->get_phi(); 02137 const std::vector<std::vector<RealGradient> >* dphi = NULL; 02138 if (cont == C_ONE) 02139 dphi = &(edge_fe->get_dphi()); 02140 02141 for (unsigned int e=0; e != elem->n_edges(); ++e) 02142 { 02143 context.edge = e; 02144 context.edge_fe_reinit(); 02145 02146 FEInterface::dofs_on_edge(elem, dim, fe_type, e, 02147 side_dofs); 02148 02149 // Some edge dofs are on nodes and already 02150 // fixed, others are free to calculate 02151 unsigned int free_dofs = 0; 02152 for (unsigned int i=0; i != side_dofs.size(); ++i) 02153 if (!dof_is_fixed[side_dofs[i]]) 02154 free_dof[free_dofs++] = i; 02155 02156 // There may be nothing to project 02157 if (!free_dofs) 02158 continue; 02159 02160 Ke.resize (free_dofs, free_dofs); Ke.zero(); 02161 Fe.resize (free_dofs); Fe.zero(); 02162 // The new edge coefficients 02163 DenseVector<Number> Uedge(free_dofs); 02164 02165 // Loop over the quadrature points 02166 for (unsigned int qp=0; qp<n_qp; qp++) 02167 { 02168 // solution at the quadrature point 02169 Number fineval = f->component(context, 02170 var_component, 02171 xyz_values[qp], 02172 system.time); 02173 // solution grad at the quadrature point 02174 Gradient finegrad; 02175 if (cont == C_ONE) 02176 finegrad = g->component(context, 02177 var_component, 02178 xyz_values[qp], 02179 system.time); 02180 02181 // Form edge projection matrix 02182 for (unsigned int sidei=0, freei=0; 02183 sidei != side_dofs.size(); ++sidei) 02184 { 02185 unsigned int i = side_dofs[sidei]; 02186 // fixed DoFs aren't test functions 02187 if (dof_is_fixed[i]) 02188 continue; 02189 for (unsigned int sidej=0, freej=0; 02190 sidej != side_dofs.size(); ++sidej) 02191 { 02192 unsigned int j = side_dofs[sidej]; 02193 if (dof_is_fixed[j]) 02194 Fe(freei) -= phi[i][qp] * phi[j][qp] * 02195 JxW[qp] * Ue(j); 02196 else 02197 Ke(freei,freej) += phi[i][qp] * 02198 phi[j][qp] * JxW[qp]; 02199 if (cont == C_ONE) 02200 { 02201 if (dof_is_fixed[j]) 02202 Fe(freei) -= ( (*dphi)[i][qp] * 02203 (*dphi)[j][qp] ) * 02204 JxW[qp] * Ue(j); 02205 else 02206 Ke(freei,freej) += ( (*dphi)[i][qp] * 02207 (*dphi)[j][qp] ) 02208 * JxW[qp]; 02209 } 02210 if (!dof_is_fixed[j]) 02211 freej++; 02212 } 02213 Fe(freei) += phi[i][qp] * fineval * JxW[qp]; 02214 if (cont == C_ONE) 02215 Fe(freei) += (finegrad * (*dphi)[i][qp] ) * 02216 JxW[qp]; 02217 freei++; 02218 } 02219 } 02220 02221 Ke.cholesky_solve(Fe, Uedge); 02222 02223 // Transfer new edge solutions to element 02224 for (unsigned int i=0; i != free_dofs; ++i) 02225 { 02226 Number &ui = Ue(side_dofs[free_dof[i]]); 02227 libmesh_assert(std::abs(ui) < TOLERANCE || 02228 std::abs(ui - Uedge(i)) < TOLERANCE); 02229 ui = Uedge(i); 02230 dof_is_fixed[side_dofs[free_dof[i]]] = true; 02231 } 02232 } 02233 } // end if (dim > 2 && cont != DISCONTINUOUS) 02234 02235 // Project any side values (edges in 2D, faces in 3D) 02236 if (dim > 1 && cont != DISCONTINUOUS) 02237 { 02238 const QBase* qsiderule = context.get_side_qrule(); 02239 const unsigned int n_qp = qsiderule->n_points(); 02240 const std::vector<Point>& xyz_values = side_fe->get_xyz(); 02241 const std::vector<Real>& JxW = side_fe->get_JxW(); 02242 02243 const std::vector<std::vector<Real> >& phi = side_fe->get_phi(); 02244 const std::vector<std::vector<RealGradient> >* dphi = NULL; 02245 if (cont == C_ONE) 02246 dphi = &(side_fe->get_dphi()); 02247 02248 for (unsigned int s=0; s != elem->n_sides(); ++s) 02249 { 02250 FEInterface::dofs_on_side(elem, dim, fe_type, s, 02251 side_dofs); 02252 02253 // Some side dofs are on nodes/edges and already 02254 // fixed, others are free to calculate 02255 unsigned int free_dofs = 0; 02256 for (unsigned int i=0; i != side_dofs.size(); ++i) 02257 if (!dof_is_fixed[side_dofs[i]]) 02258 free_dof[free_dofs++] = i; 02259 02260 // There may be nothing to project 02261 if (!free_dofs) 02262 continue; 02263 02264 Ke.resize (free_dofs, free_dofs); Ke.zero(); 02265 Fe.resize (free_dofs); Fe.zero(); 02266 // The new side coefficients 02267 DenseVector<Number> Uside(free_dofs); 02268 02269 context.side = s; 02270 context.side_fe_reinit(); 02271 02272 // Loop over the quadrature points 02273 for (unsigned int qp=0; qp<n_qp; qp++) 02274 { 02275 // solution at the quadrature point 02276 Number fineval = f->component(context, 02277 var_component, 02278 xyz_values[qp], 02279 system.time); 02280 // solution grad at the quadrature point 02281 Gradient finegrad; 02282 if (cont == C_ONE) 02283 finegrad = g->component(context, 02284 var_component, 02285 xyz_values[qp], 02286 system.time); 02287 02288 // Form side projection matrix 02289 for (unsigned int sidei=0, freei=0; 02290 sidei != side_dofs.size(); ++sidei) 02291 { 02292 unsigned int i = side_dofs[sidei]; 02293 // fixed DoFs aren't test functions 02294 if (dof_is_fixed[i]) 02295 continue; 02296 for (unsigned int sidej=0, freej=0; 02297 sidej != side_dofs.size(); ++sidej) 02298 { 02299 unsigned int j = side_dofs[sidej]; 02300 if (dof_is_fixed[j]) 02301 Fe(freei) -= phi[i][qp] * phi[j][qp] * 02302 JxW[qp] * Ue(j); 02303 else 02304 Ke(freei,freej) += phi[i][qp] * 02305 phi[j][qp] * JxW[qp]; 02306 if (cont == C_ONE) 02307 { 02308 if (dof_is_fixed[j]) 02309 Fe(freei) -= ( (*dphi)[i][qp] * 02310 (*dphi)[j][qp] ) * 02311 JxW[qp] * Ue(j); 02312 else 02313 Ke(freei,freej) += ( (*dphi)[i][qp] * 02314 (*dphi)[j][qp] ) 02315 * JxW[qp]; 02316 } 02317 if (!dof_is_fixed[j]) 02318 freej++; 02319 } 02320 Fe(freei) += (fineval * phi[i][qp]) * JxW[qp]; 02321 if (cont == C_ONE) 02322 Fe(freei) += (finegrad * (*dphi)[i][qp] ) * 02323 JxW[qp]; 02324 freei++; 02325 } 02326 } 02327 02328 Ke.cholesky_solve(Fe, Uside); 02329 02330 // Transfer new side solutions to element 02331 for (unsigned int i=0; i != free_dofs; ++i) 02332 { 02333 Number &ui = Ue(side_dofs[free_dof[i]]); 02334 libmesh_assert(std::abs(ui) < TOLERANCE || 02335 std::abs(ui - Uside(i)) < TOLERANCE); 02336 ui = Uside(i); 02337 dof_is_fixed[side_dofs[free_dof[i]]] = true; 02338 } 02339 } 02340 }// end if (dim > 1 && cont != DISCONTINUOUS) 02341 02342 // Project the interior values, finally 02343 02344 // Some interior dofs are on nodes/edges/sides and 02345 // already fixed, others are free to calculate 02346 unsigned int free_dofs = 0; 02347 for (unsigned int i=0; i != n_dofs; ++i) 02348 if (!dof_is_fixed[i]) 02349 free_dof[free_dofs++] = i; 02350 02351 // There may be nothing to project 02352 if (free_dofs) 02353 { 02354 context.elem_fe_reinit(); 02355 02356 const QBase* qrule = context.get_element_qrule(); 02357 const unsigned int n_qp = qrule->n_points(); 02358 const std::vector<Point>& xyz_values = fe->get_xyz(); 02359 const std::vector<Real>& JxW = fe->get_JxW(); 02360 02361 const std::vector<std::vector<Real> >& phi = fe->get_phi(); 02362 const std::vector<std::vector<RealGradient> >* dphi = NULL; 02363 if (cont == C_ONE) 02364 dphi = &(side_fe->get_dphi()); 02365 02366 Ke.resize (free_dofs, free_dofs); Ke.zero(); 02367 Fe.resize (free_dofs); Fe.zero(); 02368 // The new interior coefficients 02369 DenseVector<Number> Uint(free_dofs); 02370 02371 // Loop over the quadrature points 02372 for (unsigned int qp=0; qp<n_qp; qp++) 02373 { 02374 // solution at the quadrature point 02375 Number fineval = f->component(context, 02376 var_component, 02377 xyz_values[qp], 02378 system.time); 02379 // solution grad at the quadrature point 02380 Gradient finegrad; 02381 if (cont == C_ONE) 02382 finegrad = g->component(context, 02383 var_component, 02384 xyz_values[qp], 02385 system.time); 02386 02387 // Form interior projection matrix 02388 for (unsigned int i=0, freei=0; i != n_dofs; ++i) 02389 { 02390 // fixed DoFs aren't test functions 02391 if (dof_is_fixed[i]) 02392 continue; 02393 for (unsigned int j=0, freej=0; j != n_dofs; ++j) 02394 { 02395 if (dof_is_fixed[j]) 02396 Fe(freei) -= phi[i][qp] * phi[j][qp] * JxW[qp] 02397 * Ue(j); 02398 else 02399 Ke(freei,freej) += phi[i][qp] * phi[j][qp] * 02400 JxW[qp]; 02401 if (cont == C_ONE) 02402 { 02403 if (dof_is_fixed[j]) 02404 Fe(freei) -= ( (*dphi)[i][qp] * 02405 (*dphi)[j][qp] ) * JxW[qp] * 02406 Ue(j); 02407 else 02408 Ke(freei,freej) += ( (*dphi)[i][qp] * 02409 (*dphi)[j][qp] ) * 02410 JxW[qp]; 02411 } 02412 if (!dof_is_fixed[j]) 02413 freej++; 02414 } 02415 Fe(freei) += phi[i][qp] * fineval * JxW[qp]; 02416 if (cont == C_ONE) 02417 Fe(freei) += (finegrad * (*dphi)[i][qp] ) * JxW[qp]; 02418 freei++; 02419 } 02420 } 02421 Ke.cholesky_solve(Fe, Uint); 02422 02423 // Transfer new interior solutions to element 02424 for (unsigned int i=0; i != free_dofs; ++i) 02425 { 02426 Number &ui = Ue(free_dof[i]); 02427 libmesh_assert(std::abs(ui) < TOLERANCE || 02428 std::abs(ui - Uint(i)) < TOLERANCE); 02429 ui = Uint(i); 02430 dof_is_fixed[free_dof[i]] = true; 02431 } 02432 02433 } // if there are free interior dofs 02434 02435 // Make sure every DoF got reached! 02436 for (unsigned int i=0; i != n_dofs; ++i) 02437 libmesh_assert(dof_is_fixed[i]); 02438 02439 const numeric_index_type 02440 first = new_vector.first_local_index(), 02441 last = new_vector.last_local_index(); 02442 02443 // Lock the new_vector since it is shared among threads. 02444 { 02445 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 02446 02447 for (unsigned int i = 0; i < n_dofs; i++) 02448 // We may be projecting a new zero value onto 02449 // an old nonzero approximation - RHS 02450 // if (Ue(i) != 0.) 02451 if ((dof_indices[i] >= first) && 02452 (dof_indices[i] < last)) 02453 { 02454 new_vector.set(dof_indices[i], Ue(i)); 02455 } 02456 } 02457 } // end variables loop 02458 } // end elem loop 02459 } 02460 02461 02462 02463 void BoundaryProjectSolution::operator()(const ConstElemRange &range) const 02464 { 02465 // We need data to project 02466 libmesh_assert(f.get()); 02467 02475 // The dimensionality of the current mesh 02476 const unsigned int dim = system.get_mesh().mesh_dimension(); 02477 02478 // The DofMap for this system 02479 const DofMap& dof_map = system.get_dof_map(); 02480 02481 // Boundary info for the current mesh 02482 const BoundaryInfo& boundary_info = *system.get_mesh().boundary_info; 02483 02484 // The element matrix and RHS for projections. 02485 // Note that Ke is always real-valued, whereas 02486 // Fe may be complex valued if complex number 02487 // support is enabled 02488 DenseMatrix<Real> Ke; 02489 DenseVector<Number> Fe; 02490 // The new element coefficients 02491 DenseVector<Number> Ue; 02492 02493 02494 // Loop over all the variables we've been requested to project 02495 for (unsigned int v=0; v!=variables.size(); v++) 02496 { 02497 const unsigned int var = variables[v]; 02498 02499 const Variable& variable = dof_map.variable(var); 02500 02501 const FEType& fe_type = variable.type(); 02502 02503 if (fe_type.family == SCALAR) 02504 continue; 02505 02506 const unsigned int var_component = 02507 system.variable_scalar_number(var, 0); 02508 02509 // Get FE objects of the appropriate type 02510 AutoPtr<FEBase> fe (FEBase::build(dim, fe_type)); 02511 02512 // Prepare variables for projection 02513 AutoPtr<QBase> qedgerule (fe_type.default_quadrature_rule(1)); 02514 AutoPtr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1)); 02515 02516 // The values of the shape functions at the quadrature 02517 // points 02518 const std::vector<std::vector<Real> >& phi = fe->get_phi(); 02519 02520 // The gradients of the shape functions at the quadrature 02521 // points on the child element. 02522 const std::vector<std::vector<RealGradient> > *dphi = NULL; 02523 02524 const FEContinuity cont = fe->get_continuity(); 02525 02526 if (cont == C_ONE) 02527 { 02528 // We'll need gradient data for a C1 projection 02529 libmesh_assert(g.get()); 02530 02531 const std::vector<std::vector<RealGradient> >& 02532 ref_dphi = fe->get_dphi(); 02533 dphi = &ref_dphi; 02534 } 02535 02536 // The Jacobian * quadrature weight at the quadrature points 02537 const std::vector<Real>& JxW = 02538 fe->get_JxW(); 02539 02540 // The XYZ locations of the quadrature points 02541 const std::vector<Point>& xyz_values = 02542 fe->get_xyz(); 02543 02544 // The global DOF indices 02545 std::vector<dof_id_type> dof_indices; 02546 // Side/edge DOF indices 02547 std::vector<unsigned int> side_dofs; 02548 02549 // Iterate over all the elements in the range 02550 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it) 02551 { 02552 const Elem* elem = *elem_it; 02553 02554 // Per-subdomain variables don't need to be projected on 02555 // elements where they're not active 02556 if (!variable.active_on_subdomain(elem->subdomain_id())) 02557 continue; 02558 02559 // Find out which nodes, edges and sides are on a requested 02560 // boundary: 02561 std::vector<bool> is_boundary_node(elem->n_nodes(), false), 02562 is_boundary_edge(elem->n_edges(), false), 02563 is_boundary_side(elem->n_sides(), false); 02564 for (unsigned char s=0; s != elem->n_sides(); ++s) 02565 { 02566 // First see if this side has been requested 02567 const std::vector<boundary_id_type>& bc_ids = 02568 boundary_info.boundary_ids (elem, s); 02569 bool do_this_side = false; 02570 for (unsigned int i=0; i != bc_ids.size(); ++i) 02571 if (b.count(bc_ids[i])) 02572 { 02573 do_this_side = true; 02574 break; 02575 } 02576 if (!do_this_side) 02577 continue; 02578 02579 is_boundary_side[s] = true; 02580 02581 // Then see what nodes and what edges are on it 02582 for (unsigned int n=0; n != elem->n_nodes(); ++n) 02583 if (elem->is_node_on_side(n,s)) 02584 is_boundary_node[n] = true; 02585 for (unsigned int e=0; e != elem->n_edges(); ++e) 02586 if (elem->is_edge_on_side(e,s)) 02587 is_boundary_edge[e] = true; 02588 } 02589 02590 // Update the DOF indices for this element based on 02591 // the current mesh 02592 dof_map.dof_indices (elem, dof_indices, var); 02593 02594 // The number of DOFs on the element 02595 const unsigned int n_dofs = 02596 libmesh_cast_int<unsigned int>(dof_indices.size()); 02597 02598 // Fixed vs. free DoFs on edge/face projections 02599 std::vector<char> dof_is_fixed(n_dofs, false); // bools 02600 std::vector<int> free_dof(n_dofs, 0); 02601 02602 // The element type 02603 const ElemType elem_type = elem->type(); 02604 02605 // The number of nodes on the new element 02606 const unsigned int n_nodes = elem->n_nodes(); 02607 02608 // Zero the interpolated values 02609 Ue.resize (n_dofs); Ue.zero(); 02610 02611 // In general, we need a series of 02612 // projections to ensure a unique and continuous 02613 // solution. We start by interpolating boundary nodes, then 02614 // hold those fixed and project boundary edges, then hold 02615 // those fixed and project boundary faces, 02616 02617 // Interpolate node values first 02618 unsigned int current_dof = 0; 02619 for (unsigned int n=0; n!= n_nodes; ++n) 02620 { 02621 // FIXME: this should go through the DofMap, 02622 // not duplicate dof_indices code badly! 02623 const unsigned int nc = 02624 FEInterface::n_dofs_at_node (dim, fe_type, elem_type, 02625 n); 02626 if (!elem->is_vertex(n) || !is_boundary_node[n]) 02627 { 02628 current_dof += nc; 02629 continue; 02630 } 02631 if (cont == DISCONTINUOUS) 02632 { 02633 libmesh_assert_equal_to (nc, 0); 02634 } 02635 // Assume that C_ZERO elements have a single nodal 02636 // value shape function 02637 else if (cont == C_ZERO) 02638 { 02639 libmesh_assert_equal_to (nc, 1); 02640 Ue(current_dof) = f->component(var_component, 02641 elem->point(n), 02642 system.time); 02643 dof_is_fixed[current_dof] = true; 02644 current_dof++; 02645 } 02646 // The hermite element vertex shape functions are weird 02647 else if (fe_type.family == HERMITE) 02648 { 02649 Ue(current_dof) = f->component(var_component, 02650 elem->point(n), 02651 system.time); 02652 dof_is_fixed[current_dof] = true; 02653 current_dof++; 02654 Gradient grad = g->component(var_component, 02655 elem->point(n), 02656 system.time); 02657 // x derivative 02658 Ue(current_dof) = grad(0); 02659 dof_is_fixed[current_dof] = true; 02660 current_dof++; 02661 if (dim > 1) 02662 { 02663 // We'll finite difference mixed derivatives 02664 Point nxminus = elem->point(n), 02665 nxplus = elem->point(n); 02666 nxminus(0) -= TOLERANCE; 02667 nxplus(0) += TOLERANCE; 02668 Gradient gxminus = g->component(var_component, 02669 nxminus, 02670 system.time); 02671 Gradient gxplus = g->component(var_component, 02672 nxplus, 02673 system.time); 02674 // y derivative 02675 Ue(current_dof) = grad(1); 02676 dof_is_fixed[current_dof] = true; 02677 current_dof++; 02678 // xy derivative 02679 Ue(current_dof) = (gxplus(1) - gxminus(1)) 02680 / 2. / TOLERANCE; 02681 dof_is_fixed[current_dof] = true; 02682 current_dof++; 02683 02684 if (dim > 2) 02685 { 02686 // z derivative 02687 Ue(current_dof) = grad(2); 02688 dof_is_fixed[current_dof] = true; 02689 current_dof++; 02690 // xz derivative 02691 Ue(current_dof) = (gxplus(2) - gxminus(2)) 02692 / 2. / TOLERANCE; 02693 dof_is_fixed[current_dof] = true; 02694 current_dof++; 02695 // We need new points for yz 02696 Point nyminus = elem->point(n), 02697 nyplus = elem->point(n); 02698 nyminus(1) -= TOLERANCE; 02699 nyplus(1) += TOLERANCE; 02700 Gradient gyminus = g->component(var_component, 02701 nyminus, 02702 system.time); 02703 Gradient gyplus = g->component(var_component, 02704 nyplus, 02705 system.time); 02706 // xz derivative 02707 Ue(current_dof) = (gyplus(2) - gyminus(2)) 02708 / 2. / TOLERANCE; 02709 dof_is_fixed[current_dof] = true; 02710 current_dof++; 02711 // Getting a 2nd order xyz is more tedious 02712 Point nxmym = elem->point(n), 02713 nxmyp = elem->point(n), 02714 nxpym = elem->point(n), 02715 nxpyp = elem->point(n); 02716 nxmym(0) -= TOLERANCE; 02717 nxmym(1) -= TOLERANCE; 02718 nxmyp(0) -= TOLERANCE; 02719 nxmyp(1) += TOLERANCE; 02720 nxpym(0) += TOLERANCE; 02721 nxpym(1) -= TOLERANCE; 02722 nxpyp(0) += TOLERANCE; 02723 nxpyp(1) += TOLERANCE; 02724 Gradient gxmym = g->component(var_component, 02725 nxmym, 02726 system.time); 02727 Gradient gxmyp = g->component(var_component, 02728 nxmyp, 02729 system.time); 02730 Gradient gxpym = g->component(var_component, 02731 nxpym, 02732 system.time); 02733 Gradient gxpyp = g->component(var_component, 02734 nxpyp, 02735 system.time); 02736 Number gxzplus = (gxpyp(2) - gxmyp(2)) 02737 / 2. / TOLERANCE; 02738 Number gxzminus = (gxpym(2) - gxmym(2)) 02739 / 2. / TOLERANCE; 02740 // xyz derivative 02741 Ue(current_dof) = (gxzplus - gxzminus) 02742 / 2. / TOLERANCE; 02743 dof_is_fixed[current_dof] = true; 02744 current_dof++; 02745 } 02746 } 02747 } 02748 // Assume that other C_ONE elements have a single nodal 02749 // value shape function and nodal gradient component 02750 // shape functions 02751 else if (cont == C_ONE) 02752 { 02753 libmesh_assert_equal_to (nc, 1 + dim); 02754 Ue(current_dof) = f->component(var_component, 02755 elem->point(n), 02756 system.time); 02757 dof_is_fixed[current_dof] = true; 02758 current_dof++; 02759 Gradient grad = g->component(var_component, 02760 elem->point(n), 02761 system.time); 02762 for (unsigned int i=0; i!= dim; ++i) 02763 { 02764 Ue(current_dof) = grad(i); 02765 dof_is_fixed[current_dof] = true; 02766 current_dof++; 02767 } 02768 } 02769 else 02770 libmesh_error(); 02771 } 02772 02773 // In 3D, project any edge values next 02774 if (dim > 2 && cont != DISCONTINUOUS) 02775 for (unsigned int e=0; e != elem->n_edges(); ++e) 02776 { 02777 if (!is_boundary_edge[e]) 02778 continue; 02779 02780 FEInterface::dofs_on_edge(elem, dim, fe_type, e, 02781 side_dofs); 02782 02783 // Some edge dofs are on nodes and already 02784 // fixed, others are free to calculate 02785 unsigned int free_dofs = 0; 02786 for (unsigned int i=0; i != side_dofs.size(); ++i) 02787 if (!dof_is_fixed[side_dofs[i]]) 02788 free_dof[free_dofs++] = i; 02789 02790 // There may be nothing to project 02791 if (!free_dofs) 02792 continue; 02793 02794 Ke.resize (free_dofs, free_dofs); Ke.zero(); 02795 Fe.resize (free_dofs); Fe.zero(); 02796 // The new edge coefficients 02797 DenseVector<Number> Uedge(free_dofs); 02798 02799 // Initialize FE data on the edge 02800 fe->attach_quadrature_rule (qedgerule.get()); 02801 fe->edge_reinit (elem, e); 02802 const unsigned int n_qp = qedgerule->n_points(); 02803 02804 // Loop over the quadrature points 02805 for (unsigned int qp=0; qp<n_qp; qp++) 02806 { 02807 // solution at the quadrature point 02808 Number fineval = f->component(var_component, 02809 xyz_values[qp], 02810 system.time); 02811 // solution grad at the quadrature point 02812 Gradient finegrad; 02813 if (cont == C_ONE) 02814 finegrad = g->component(var_component, 02815 xyz_values[qp], 02816 system.time); 02817 02818 // Form edge projection matrix 02819 for (unsigned int sidei=0, freei=0; 02820 sidei != side_dofs.size(); ++sidei) 02821 { 02822 unsigned int i = side_dofs[sidei]; 02823 // fixed DoFs aren't test functions 02824 if (dof_is_fixed[i]) 02825 continue; 02826 for (unsigned int sidej=0, freej=0; 02827 sidej != side_dofs.size(); ++sidej) 02828 { 02829 unsigned int j = side_dofs[sidej]; 02830 if (dof_is_fixed[j]) 02831 Fe(freei) -= phi[i][qp] * phi[j][qp] * 02832 JxW[qp] * Ue(j); 02833 else 02834 Ke(freei,freej) += phi[i][qp] * 02835 phi[j][qp] * JxW[qp]; 02836 if (cont == C_ONE) 02837 { 02838 if (dof_is_fixed[j]) 02839 Fe(freei) -= ((*dphi)[i][qp] * 02840 (*dphi)[j][qp]) * 02841 JxW[qp] * Ue(j); 02842 else 02843 Ke(freei,freej) += ((*dphi)[i][qp] * 02844 (*dphi)[j][qp]) 02845 * JxW[qp]; 02846 } 02847 if (!dof_is_fixed[j]) 02848 freej++; 02849 } 02850 Fe(freei) += phi[i][qp] * fineval * JxW[qp]; 02851 if (cont == C_ONE) 02852 Fe(freei) += (finegrad * (*dphi)[i][qp]) * 02853 JxW[qp]; 02854 freei++; 02855 } 02856 } 02857 02858 Ke.cholesky_solve(Fe, Uedge); 02859 02860 // Transfer new edge solutions to element 02861 for (unsigned int i=0; i != free_dofs; ++i) 02862 { 02863 Number &ui = Ue(side_dofs[free_dof[i]]); 02864 libmesh_assert(std::abs(ui) < TOLERANCE || 02865 std::abs(ui - Uedge(i)) < TOLERANCE); 02866 ui = Uedge(i); 02867 dof_is_fixed[side_dofs[free_dof[i]]] = true; 02868 } 02869 } 02870 02871 // Project any side values (edges in 2D, faces in 3D) 02872 if (dim > 1 && cont != DISCONTINUOUS) 02873 for (unsigned int s=0; s != elem->n_sides(); ++s) 02874 { 02875 if (!is_boundary_side[s]) 02876 continue; 02877 02878 FEInterface::dofs_on_side(elem, dim, fe_type, s, 02879 side_dofs); 02880 02881 // Some side dofs are on nodes/edges and already 02882 // fixed, others are free to calculate 02883 unsigned int free_dofs = 0; 02884 for (unsigned int i=0; i != side_dofs.size(); ++i) 02885 if (!dof_is_fixed[side_dofs[i]]) 02886 free_dof[free_dofs++] = i; 02887 02888 // There may be nothing to project 02889 if (!free_dofs) 02890 continue; 02891 02892 Ke.resize (free_dofs, free_dofs); Ke.zero(); 02893 Fe.resize (free_dofs); Fe.zero(); 02894 // The new side coefficients 02895 DenseVector<Number> Uside(free_dofs); 02896 02897 // Initialize FE data on the side 02898 fe->attach_quadrature_rule (qsiderule.get()); 02899 fe->reinit (elem, s); 02900 const unsigned int n_qp = qsiderule->n_points(); 02901 02902 // Loop over the quadrature points 02903 for (unsigned int qp=0; qp<n_qp; qp++) 02904 { 02905 // solution at the quadrature point 02906 Number fineval = f->component(var_component, 02907 xyz_values[qp], 02908 system.time); 02909 // solution grad at the quadrature point 02910 Gradient finegrad; 02911 if (cont == C_ONE) 02912 finegrad = g->component(var_component, 02913 xyz_values[qp], 02914 system.time); 02915 02916 // Form side projection matrix 02917 for (unsigned int sidei=0, freei=0; 02918 sidei != side_dofs.size(); ++sidei) 02919 { 02920 unsigned int i = side_dofs[sidei]; 02921 // fixed DoFs aren't test functions 02922 if (dof_is_fixed[i]) 02923 continue; 02924 for (unsigned int sidej=0, freej=0; 02925 sidej != side_dofs.size(); ++sidej) 02926 { 02927 unsigned int j = side_dofs[sidej]; 02928 if (dof_is_fixed[j]) 02929 Fe(freei) -= phi[i][qp] * phi[j][qp] * 02930 JxW[qp] * Ue(j); 02931 else 02932 Ke(freei,freej) += phi[i][qp] * 02933 phi[j][qp] * JxW[qp]; 02934 if (cont == C_ONE) 02935 { 02936 if (dof_is_fixed[j]) 02937 Fe(freei) -= ((*dphi)[i][qp] * 02938 (*dphi)[j][qp]) * 02939 JxW[qp] * Ue(j); 02940 else 02941 Ke(freei,freej) += ((*dphi)[i][qp] * 02942 (*dphi)[j][qp]) 02943 * JxW[qp]; 02944 } 02945 if (!dof_is_fixed[j]) 02946 freej++; 02947 } 02948 Fe(freei) += (fineval * phi[i][qp]) * JxW[qp]; 02949 if (cont == C_ONE) 02950 Fe(freei) += (finegrad * (*dphi)[i][qp]) * 02951 JxW[qp]; 02952 freei++; 02953 } 02954 } 02955 02956 Ke.cholesky_solve(Fe, Uside); 02957 02958 // Transfer new side solutions to element 02959 for (unsigned int i=0; i != free_dofs; ++i) 02960 { 02961 Number &ui = Ue(side_dofs[free_dof[i]]); 02962 libmesh_assert(std::abs(ui) < TOLERANCE || 02963 std::abs(ui - Uside(i)) < TOLERANCE); 02964 ui = Uside(i); 02965 dof_is_fixed[side_dofs[free_dof[i]]] = true; 02966 } 02967 } 02968 02969 const dof_id_type 02970 first = new_vector.first_local_index(), 02971 last = new_vector.last_local_index(); 02972 02973 // Lock the new_vector since it is shared among threads. 02974 { 02975 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 02976 02977 for (unsigned int i = 0; i < n_dofs; i++) 02978 if (dof_is_fixed[i] && 02979 (dof_indices[i] >= first) && 02980 (dof_indices[i] < last)) 02981 { 02982 new_vector.set(dof_indices[i], Ue(i)); 02983 } 02984 } 02985 } // end elem loop 02986 } // end variables loop 02987 } 02988 02989 02990 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: