system_projection.C

Go to the documentation of this file.
00001 // $Id: system_projection.C 3530 2009-10-30 23:03:08Z roystgnr $
00002 
00003 // The libMesh Finite Element Library.
00004 // Copyright (C) 2002-2008 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00005   
00006 // This library is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 2.1 of the License, or (at your option) any later version.
00010   
00011 // This library is distributed in the hope that it will be useful,
00012 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00014 // Lesser General Public License for more details.
00015   
00016 // You should have received a copy of the GNU Lesser General Public
00017 // License along with this library; if not, write to the Free Software
00018 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019 
00020 
00021 
00022 // C++ includes
00023 #include <vector>
00024 
00025 // Local includes
00026 #include "system.h"
00027 #include "mesh.h"
00028 #include "elem.h"
00029 #include "dof_map.h"
00030 #include "fe_interface.h"
00031 #include "numeric_vector.h"
00032 #include "libmesh_logging.h"
00033 
00034 #include "dense_matrix.h"
00035 #include "fe_base.h"
00036 #include "quadrature_gauss.h"
00037 #include "dense_vector.h"
00038 #include "threads.h"
00039 
00040 
00041 // ------------------------------------------------------------
00042 // System implementation
00043 void System::project_vector (NumericVector<Number>& vector) const
00044 {
00045   // Create a copy of the vector, which currently
00046   // contains the old data.
00047   AutoPtr<NumericVector<Number> >
00048     old_vector (vector.clone());
00049 
00050   // Project the old vector to the new vector
00051   this->project_vector (*old_vector, vector);
00052 }
00053 
00054 
00060 void System::project_vector (const NumericVector<Number>& old_v,
00061                              NumericVector<Number>& new_v) const
00062 {
00063   START_LOG ("project_vector()", "System");
00064 
00071   new_v.clear();
00072 
00073 #ifdef LIBMESH_ENABLE_AMR
00074 
00075   // Resize the new vector and get a serial version.
00076   NumericVector<Number> *new_vector_ptr = NULL;
00077   AutoPtr<NumericVector<Number> > new_vector_built;
00078   NumericVector<Number> *local_old_vector;
00079   AutoPtr<NumericVector<Number> > local_old_vector_built;
00080   const NumericVector<Number> *old_vector_ptr = NULL;
00081 
00082   ConstElemRange active_local_elem_range 
00083     (this->get_mesh().active_local_elements_begin(),
00084      this->get_mesh().active_local_elements_end());
00085 
00086   // If the old vector was uniprocessor, make the new
00087   // vector uniprocessor
00088   if (old_v.type() == SERIAL)
00089     {
00090       new_v.init (this->n_dofs(), false, SERIAL);
00091       new_vector_ptr = &new_v;
00092       old_vector_ptr = &old_v;
00093     }
00094 
00095   // Otherwise it is a parallel, distributed vector, which
00096   // we need to localize.
00097   else if (old_v.type() == PARALLEL)
00098     {
00099       // Build a send list for efficient localization
00100       BuildProjectionList projection_list(*this);
00101       Threads::parallel_reduce (active_local_elem_range, 
00102                                 projection_list);
00103 
00104       // Create a sorted, unique send_list
00105       projection_list.unique();
00106 
00107       new_v.init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
00108       new_vector_built = NumericVector<Number>::build();
00109       local_old_vector_built = NumericVector<Number>::build();
00110       new_vector_ptr = new_vector_built.get();
00111       local_old_vector = local_old_vector_built.get();
00112       new_vector_ptr->init(this->n_dofs(), false, SERIAL);
00113       local_old_vector->init(old_v.size(), false, SERIAL);
00114       old_v.localize(*local_old_vector, projection_list.send_list);
00115       local_old_vector->close();
00116       old_vector_ptr = local_old_vector;
00117     }
00118   else if (old_v.type() == GHOSTED)
00119     {
00120       // Build a send list for efficient localization
00121       BuildProjectionList projection_list(*this);
00122       Threads::parallel_reduce (active_local_elem_range, 
00123                                 projection_list);
00124 
00125       // Create a sorted, unique send_list
00126       projection_list.unique();
00127 
00128       new_v.init (this->n_dofs(), this->n_local_dofs(),
00129                   this->get_dof_map().get_send_list(), false, GHOSTED);
00130 
00131       local_old_vector_built = NumericVector<Number>::build();
00132       new_vector_ptr = &new_v;
00133       local_old_vector = local_old_vector_built.get();
00134       local_old_vector->init(old_v.size(), old_v.local_size(),
00135                              projection_list.send_list, false, GHOSTED);
00136       old_v.localize(*local_old_vector, projection_list.send_list);
00137       local_old_vector->close();
00138       old_vector_ptr = local_old_vector;
00139     }
00140   else // unknown old_v.type()
00141     {
00142       std::cerr << "ERROR: Unknown old_v.type() == " << old_v.type() 
00143                 << std::endl;
00144       libmesh_error();
00145     }
00146   
00147   // Note that the above will have zeroed the new_vector.
00148   // Just to be sure, assert that new_vector_ptr and old_vector_ptr
00149   // were successfully set before trying to deref them.
00150   libmesh_assert(new_vector_ptr);
00151   libmesh_assert(old_vector_ptr);
00152 
00153   NumericVector<Number> &new_vector = *new_vector_ptr;
00154   const NumericVector<Number> &old_vector = *old_vector_ptr;
00155     
00156   Threads::parallel_for (active_local_elem_range,
00157                          ProjectVector(*this,
00158                                        old_vector,
00159                                        new_vector)
00160                          );
00161 
00162   // Copy the SCALAR dofs from old_vector to new_vector
00163   // Note: We assume that all SCALAR dofs are on the
00164   // processor with highest ID
00165   if(libMesh::processor_id() == (libMesh::n_processors()-1))
00166   {
00167     const DofMap& dof_map = this->get_dof_map();
00168     for (unsigned int var=0; var<this->n_vars(); var++)
00169       if(this->variable(var).type().family == SCALAR)
00170         {
00171           // We can just map SCALAR dofs directly across
00172           std::vector<unsigned int> new_SCALAR_indices, old_SCALAR_indices;
00173           dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false);
00174           dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true);
00175           const unsigned int new_n_dofs = new_SCALAR_indices.size();
00176 
00177           for (unsigned int i=0; i<new_n_dofs; i++)
00178           {
00179             new_vector.set( new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]) );
00180           }
00181         }
00182   }
00183 
00184   new_vector.close();
00185 
00186   // If the old vector was serial, we probably need to send our values
00187   // to other processors
00188   //
00189   // FIXME: I'm not sure how to make a NumericVector do that without
00190   // creating a temporary parallel vector to use localize! - RHS
00191   if (old_v.type() == SERIAL)
00192     {
00193       AutoPtr<NumericVector<Number> > dist_v = NumericVector<Number>::build();
00194       dist_v->init(this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
00195       dist_v->close();
00196     
00197       for (unsigned int i=0; i!=dist_v->size(); i++)
00198         if (new_vector(i) != 0.0)
00199           dist_v->set(i, new_vector(i));
00200 
00201       dist_v->close();
00202 
00203       dist_v->localize (new_v, this->get_dof_map().get_send_list());
00204       new_v.close();
00205     }
00206   // If the old vector was parallel, we need to update it
00207   // and free the localized copies
00208   else if (old_v.type() == PARALLEL)
00209     {
00210       // We may have to set dof values that this processor doesn't
00211       // own in certain special cases, like LAGRANGE FIRST or
00212       // HERMITE THIRD elements on second-order meshes
00213       for (unsigned int i=0; i!=new_v.size(); i++)
00214         if (new_vector(i) != 0.0)
00215           new_v.set(i, new_vector(i));
00216       new_v.close();
00217     }
00218 
00219   this->get_dof_map().enforce_constraints_exactly(*this, &new_v);
00220 
00221 #else
00222 
00223   // AMR is disabled: simply copy the vector
00224   new_v = old_v;
00225   
00226 #endif // #ifdef LIBMESH_ENABLE_AMR
00227 
00228   STOP_LOG("project_vector()", "System");
00229 }
00230 
00231 
00232 
00237 void System::project_solution (Number fptr(const Point& p,
00238                                            const Parameters& parameters,
00239                                            const std::string& sys_name,
00240                                            const std::string& unknown_name),
00241                                Gradient gptr(const Point& p,
00242                                              const Parameters& parameters,
00243                                              const std::string& sys_name,
00244                                              const std::string& unknown_name),
00245                                Parameters& parameters) const
00246 {
00247   this->project_vector(fptr, gptr, parameters, *solution);
00248 
00249   solution->localize(*current_local_solution);
00250 }
00251 
00252 
00253 
00258 void System::project_vector (Number fptr(const Point& p,
00259                                          const Parameters& parameters,
00260                                          const std::string& sys_name,
00261                                          const std::string& unknown_name),
00262                              Gradient gptr(const Point& p,
00263                                            const Parameters& parameters,
00264                                            const std::string& sys_name,
00265                                            const std::string& unknown_name),
00266                              Parameters& parameters,
00267                              NumericVector<Number>& new_vector) const
00268 {
00269   START_LOG ("project_vector()", "System");
00270 
00271   Threads::parallel_for (ConstElemRange (this->get_mesh().active_local_elements_begin(),
00272                                          this->get_mesh().active_local_elements_end(),
00273                                          1000),
00274                          ProjectSolution(*this,
00275                                          fptr,
00276                                          gptr,
00277                                          parameters,
00278                                          new_vector)
00279                          );
00280 
00281   // Also, load values into the SCALAR dofs
00282   // Note: We assume that all SCALAR dofs are on the
00283   // processor with highest ID
00284   if(libMesh::processor_id() == (libMesh::n_processors()-1))
00285   {
00286     const DofMap& dof_map = this->get_dof_map();
00287     for (unsigned int var=0; var<this->n_vars(); var++)
00288       if(this->variable(var).type().family == SCALAR)
00289         {
00290           std::vector<unsigned int> SCALAR_indices;
00291           dof_map.SCALAR_dof_indices (SCALAR_indices, var);
00292           const unsigned int n_SCALAR_dofs = SCALAR_indices.size();
00293 
00294           for (unsigned int i=0; i<n_SCALAR_dofs; i++)
00295           {
00296             const unsigned int index = SCALAR_indices[i];
00297 
00298             // We pass the point (i,0,0) to the fptr to distinguish
00299             // the different scalars within the SCALAR variable
00300             Point p_i(i,0,0);
00301             new_vector.set( index, fptr(p_i,
00302                                         parameters,
00303                                         this->name(),
00304                                         this->variable_name(var))
00305                           );
00306           }
00307         }
00308   }
00309 
00310   new_vector.close();
00311 
00312 #ifdef LIBMESH_ENABLE_AMR
00313   this->get_dof_map().enforce_constraints_exactly(*this, &new_vector);
00314 #endif
00315 
00316   STOP_LOG("project_vector()", "System");
00317 }
00318 
00319 
00320 #ifndef LIBMESH_ENABLE_AMR
00321 void System::ProjectVector::operator()(const ConstElemRange &) const
00322 {
00323   libmesh_error();
00324 #else
00325 void System::ProjectVector::operator()(const ConstElemRange &range) const
00326 {
00327   START_LOG ("operator()","ProjectVector");
00328 
00329   // A vector for Lagrange element interpolation, indicating if we
00330   // have visited a DOF yet.  Note that this is thread-local storage,
00331   // hence shared DOFS that live on thread boundaries may be doubly
00332   // computed.  It is expected that this will be more efficient
00333   // than locking a thread-global version of already_done, though.
00334   //
00335   // FIXME: we should use this for non-Lagrange coarsening, too
00336   std::vector<bool> already_done (system.n_dofs(), false);
00337 
00338   // The number of variables in this system
00339   const unsigned int n_variables = system.n_vars();
00340 
00341   // The dimensionality of the current mesh
00342   const unsigned int dim = system.get_mesh().mesh_dimension();
00343 
00344   // The DofMap for this system
00345   const DofMap& dof_map = system.get_dof_map();
00346 
00347   // The element matrix and RHS for projections.
00348   // Note that Ke is always real-valued, whereas
00349   // Fe may be complex valued if complex number
00350   // support is enabled
00351   DenseMatrix<Real> Ke;
00352   DenseVector<Number> Fe;
00353   // The new element coefficients
00354   DenseVector<Number> Ue;
00355 
00356 
00357   // Loop over all the variables in the system
00358   for (unsigned int var=0; var<n_variables; var++)
00359     {
00360       if (dof_map.variable(var).type().family == SCALAR)
00361         continue;
00362 
00363       // Get FE objects of the appropriate type
00364       const FEType& base_fe_type = dof_map.variable_type(var);     
00365       AutoPtr<FEBase> fe (FEBase::build(dim, base_fe_type));      
00366       AutoPtr<FEBase> fe_coarse (FEBase::build(dim, base_fe_type));      
00367 
00368       // Create FE objects with potentially different p_level
00369       FEType fe_type, temp_fe_type;
00370 
00371       // Prepare variables for non-Lagrange projection
00372       AutoPtr<QBase> qrule     (base_fe_type.default_quadrature_rule(dim));
00373       AutoPtr<QBase> qedgerule (base_fe_type.default_quadrature_rule(1));
00374       AutoPtr<QBase> qsiderule (base_fe_type.default_quadrature_rule(dim-1));
00375       std::vector<Point> coarse_qpoints;
00376 
00377       // The values of the shape functions at the quadrature
00378       // points
00379       const std::vector<std::vector<Real> >& phi_values =
00380         fe->get_phi();
00381       const std::vector<std::vector<Real> >& phi_coarse =
00382         fe_coarse->get_phi();
00383 
00384       // The gradients of the shape functions at the quadrature
00385       // points on the child element.
00386       const std::vector<std::vector<RealGradient> > *dphi_values =
00387         NULL;
00388       const std::vector<std::vector<RealGradient> > *dphi_coarse =
00389         NULL;
00390 
00391       const FEContinuity cont = fe->get_continuity();
00392 
00393       if (cont == C_ONE)
00394         {
00395           const std::vector<std::vector<RealGradient> >&
00396             ref_dphi_values = fe->get_dphi();
00397           dphi_values = &ref_dphi_values;
00398           const std::vector<std::vector<RealGradient> >&
00399             ref_dphi_coarse = fe_coarse->get_dphi();
00400           dphi_coarse = &ref_dphi_coarse;
00401         }
00402 
00403       // The Jacobian * quadrature weight at the quadrature points
00404       const std::vector<Real>& JxW =
00405         fe->get_JxW();
00406      
00407       // The XYZ locations of the quadrature points on the
00408       // child element
00409       const std::vector<Point>& xyz_values =
00410         fe->get_xyz();
00411 
00412 
00413       // The global DOF indices
00414       std::vector<unsigned int> new_dof_indices, old_dof_indices;
00415       // Side/edge DOF indices
00416       std::vector<unsigned int> new_side_dofs, old_side_dofs;
00417    
00418       // Iterate over the elements in the range
00419       for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it)
00420         {
00421           const Elem* elem = *elem_it;
00422           const Elem* parent = elem->parent();
00423 
00424           // Adjust the FE type for p-refined elements
00425           fe_type = base_fe_type;
00426           fe_type.order = static_cast<Order>(fe_type.order +
00427                                              elem->p_level());
00428 
00429           // We may need to remember the parent's p_level
00430           unsigned int old_parent_level = 0;
00431 
00432           // Update the DOF indices for this element based on
00433           // the new mesh
00434           dof_map.dof_indices (elem, new_dof_indices, var);
00435 
00436           // The number of DOFs on the new element
00437           const unsigned int new_n_dofs = new_dof_indices.size();
00438 
00439           // Fixed vs. free DoFs on edge/face projections
00440           std::vector<char> dof_is_fixed(new_n_dofs, false); // bools
00441           std::vector<int> free_dof(new_n_dofs, 0);
00442 
00443           // The element type
00444           const ElemType elem_type = elem->type();
00445 
00446           // The number of nodes on the new element
00447           const unsigned int n_nodes = elem->n_nodes();
00448 
00449           // Zero the interpolated values
00450           Ue.resize (new_n_dofs); Ue.zero();
00451 
00452           // Update the DOF indices based on the old mesh.
00453           // This is done in one of three ways:
00454           // 1.) If the child was just refined then it was not
00455           //     active in the previous mesh & hence has no solution
00456           //     values on it.  In this case simply project or
00457           //     interpolate the solution from the parent, who was
00458           //     active in the previous mesh
00459           // 2.) If the child was just coarsened, obtaining a
00460           //     well-defined solution may require doing independent
00461           //     projections on nodes, edges, faces, and interiors
00462           // 3.) If the child was active in the previous
00463           //     mesh, we can just copy coefficients directly
00464           if (elem->refinement_flag() == Elem::JUST_REFINED)
00465             {
00466               libmesh_assert (parent != NULL);
00467               old_parent_level = parent->p_level();
00468 
00469               // We may have done p refinement or coarsening as well;
00470               // if so then we need to reset the parent's p level
00471               // so we can get the right DoFs from it
00472               if (elem->p_refinement_flag() == Elem::JUST_REFINED)
00473                 {
00474                   libmesh_assert(elem->p_level() > 0);
00475                   (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() - 1);
00476                 }
00477               else if (elem->p_refinement_flag() == Elem::JUST_COARSENED)
00478                 {
00479                   (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() + 1);
00480                 }
00481          
00482               dof_map.old_dof_indices (parent, old_dof_indices, var);
00483             }
00484           else
00485             {
00486               dof_map.old_dof_indices (elem, old_dof_indices, var);
00487 
00488               if (elem->p_refinement_flag() == Elem::DO_NOTHING)
00489                 libmesh_assert (old_dof_indices.size() == new_n_dofs);
00490               if (elem->p_refinement_flag() == Elem::JUST_COARSENED)
00491                 libmesh_assert (elem->has_children());
00492             }
00493 
00494           unsigned int old_n_dofs = old_dof_indices.size();
00495 
00496           if (fe_type.family != LAGRANGE) {
00497 
00498             // For refined non-Lagrange elements, we do an L2
00499             // projection
00500             // FIXME: this will be a suboptimal and ill-defined
00501             // result if we're using non-nested finite element
00502             // spaces or if we're on a p-coarsened element!
00503             if (elem->refinement_flag() == Elem::JUST_REFINED)
00504               {
00505                 // Update the fe object based on the current child
00506                 fe->attach_quadrature_rule (qrule.get());
00507                 fe->reinit (elem);
00508           
00509                 // The number of quadrature points on the child
00510                 const unsigned int n_qp = qrule->n_points();
00511 
00512                 FEInterface::inverse_map (dim, fe_type, parent,
00513                                           xyz_values, coarse_qpoints);
00514 
00515                 fe_coarse->reinit(parent, &coarse_qpoints);
00516 
00517                 // Reinitialize the element matrix and vector for
00518                 // the current element.  Note that this will zero them
00519                 // before they are summed.
00520                 Ke.resize (new_n_dofs, new_n_dofs); Ke.zero();
00521                 Fe.resize (new_n_dofs); Fe.zero();
00522           
00523           
00524                 // Loop over the quadrature points
00525                 for (unsigned int qp=0; qp<n_qp; qp++)
00526                   {
00527                     // The solution value at the quadrature point             
00528                     Number val = libMesh::zero;
00529 
00530                     // Sum the function values * the DOF values
00531                     // at the point of interest to get the function value
00532                     // (Note that the # of DOFs on the parent need not be the
00533                     //  same as on the child!)
00534                     for (unsigned int i=0; i<old_n_dofs; i++)
00535                       {
00536                         val += (old_vector(old_dof_indices[i])*
00537                                 phi_coarse[i][qp]);
00538                       }
00539 
00540                     // Now \p val contains the solution value of variable
00541                     // \p var at the qp'th quadrature point on the child
00542                     // element.  It has been interpolated from the parent
00543                     // in case the child was just refined.  Next we will
00544                     // construct the L2-projection matrix for the element.
00545 
00546                     // Construct the Mass Matrix
00547                     for (unsigned int i=0; i<new_n_dofs; i++)
00548                       for (unsigned int j=0; j<new_n_dofs; j++)
00549                         Ke(i,j) += JxW[qp]*phi_values[i][qp]*phi_values[j][qp];
00550 
00551                     // Construct the RHS
00552                     for (unsigned int i=0; i<new_n_dofs; i++)
00553                       Fe(i) += JxW[qp]*phi_values[i][qp]*val;
00554               
00555                   } // end qp loop
00556 
00557                 Ke.cholesky_solve(Fe, Ue);
00558 
00559                 // Fix up the parent's p level in case we changed it
00560                 (const_cast<Elem *>(parent))->hack_p_level(old_parent_level);
00561               }
00562             else if (elem->refinement_flag() == Elem::JUST_COARSENED)
00563               {
00564                 FEBase::coarsened_dof_values(old_vector, dof_map,
00565                                              elem, Ue, var, true);
00566               }
00567             // For unrefined uncoarsened elements, we just copy DoFs
00568             else
00569               {
00570                 // FIXME - I'm sure this function would be about half
00571                 // the size if anyone ever figures out how to improve
00572                 // the DofMap interface... - RHS
00573                 if (elem->p_refinement_flag() == Elem::JUST_REFINED)
00574                   {
00575                     libmesh_assert (elem->p_level() > 0);
00576                     // P refinement of non-hierarchic bases will
00577                     // require a whole separate code path
00578                     libmesh_assert (fe->is_hierarchic());
00579                     temp_fe_type = fe_type;
00580                     temp_fe_type.order =
00581                       static_cast<Order>(temp_fe_type.order - 1);
00582                     unsigned int old_index = 0, new_index = 0;
00583                     for (unsigned int n=0; n != n_nodes; ++n)
00584                       {
00585                         const unsigned int nc =
00586                           FEInterface::n_dofs_at_node (dim, temp_fe_type,
00587                                                        elem_type, n);
00588                         for (unsigned int i=0; i != nc; ++i)
00589                           {
00590                             Ue(new_index + i) =
00591                               old_vector(old_dof_indices[old_index++]);
00592                           }
00593                         new_index +=
00594                           FEInterface::n_dofs_at_node (dim, fe_type,
00595                                                        elem_type, n);
00596                       }
00597                   }
00598                 else if (elem->p_refinement_flag() ==
00599                          Elem::JUST_COARSENED)
00600                   {
00601                     // P coarsening of non-hierarchic bases will
00602                     // require a whole separate code path
00603                     libmesh_assert (fe->is_hierarchic());
00604                     temp_fe_type = fe_type;
00605                     temp_fe_type.order =
00606                       static_cast<Order>(temp_fe_type.order + 1);
00607                     unsigned int old_index = 0, new_index = 0;
00608                     for (unsigned int n=0; n != n_nodes; ++n)
00609                       {
00610                         const unsigned int nc =
00611                           FEInterface::n_dofs_at_node (dim, fe_type,
00612                                                        elem_type, n);
00613                         for (unsigned int i=0; i != nc; ++i)
00614                           {
00615                             Ue(new_index++) =
00616                               old_vector(old_dof_indices[old_index+i]);
00617                           }
00618                         old_index +=
00619                           FEInterface::n_dofs_at_node (dim, temp_fe_type,
00620                                                        elem_type, n);
00621                       }
00622                   }
00623                 else
00624                   // If there's no p refinement, we can copy every DoF
00625                   for (unsigned int i=0; i<new_n_dofs; i++)
00626                     Ue(i) = old_vector(old_dof_indices[i]);
00627               }
00628           } 
00629           else { // fe type is Lagrange
00630             // Loop over the DOFs on the element
00631             for (unsigned int new_local_dof=0;
00632                  new_local_dof<new_n_dofs; new_local_dof++)
00633               {
00634                 // The global DOF number for this local DOF
00635                 const unsigned int new_global_dof =
00636                   new_dof_indices[new_local_dof];
00637 
00638                 // The global DOF might lie outside of the bounds of a
00639                 // distributed vector.  Check for that and possibly continue
00640                 // the loop
00641                 if ((new_global_dof <  new_vector.first_local_index()) ||
00642                     (new_global_dof >= new_vector.last_local_index()))
00643                   continue;
00644 
00645                 // We might have already computed the solution for this DOF.
00646                 // This is likely in the case of a shared node, particularly
00647                 // at the corners of an element.  Check to see if that is the
00648                 // case
00649                 if (already_done[new_global_dof] == true)
00650                   continue;
00651 
00652                 already_done[new_global_dof] = true;
00653               
00654                 if (elem->refinement_flag() == Elem::JUST_REFINED)
00655                   {
00656                     // The location of the child's node on the parent element
00657                     const Point point =
00658                       FEInterface::inverse_map (dim, fe_type, parent,
00659                                                 elem->point(new_local_dof));
00660                   
00661                     // Sum the function values * the DOF values
00662                     // at the point of interest to get the function value
00663                     // (Note that the # of DOFs on the parent need not be the
00664                     //  same as on the child!)
00665                     for (unsigned int old_local_dof=0;
00666                          old_local_dof<old_n_dofs; old_local_dof++)
00667                       {
00668                         const unsigned int old_global_dof =
00669                           old_dof_indices[old_local_dof];
00670                       
00671                         Ue(new_local_dof) += 
00672                           (old_vector(old_global_dof)*
00673                           FEInterface::shape(dim, fe_type, parent,
00674                                              old_local_dof, point));
00675                       }
00676                   }
00677                 else
00678                   {
00679                     // Get the old global DOF index
00680                     const unsigned int old_global_dof =
00681                       old_dof_indices[new_local_dof];
00682 
00683                     Ue(new_local_dof) = old_vector(old_global_dof);
00684                   }
00685               } // end local DOF loop
00686 
00687             // We may have to clean up a parent's p_level
00688             if (elem->refinement_flag() == Elem::JUST_REFINED)
00689               (const_cast<Elem *>(parent))->hack_p_level(old_parent_level);
00690           }  // end fe_type if()
00691 
00692           // Lock the new_vector since it is shared among threads.
00693           {
00694             Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00695 
00696             for (unsigned int i = 0; i < new_n_dofs; i++) 
00697               if (Ue(i) != 0.)
00698                 new_vector.set(new_dof_indices[i], Ue(i));
00699           }
00700         }  // end elem loop
00701     } // end variables loop
00702 
00703   STOP_LOG ("operator()","ProjectVector");
00704 
00705 #endif // LIBMESH_ENABLE_AMR
00706 }
00707 
00708 
00709 
00710 void System::BuildProjectionList::unique()
00711 {
00712   // Sort the send list.  After this duplicated
00713   // elements will be adjacent in the vector
00714   std::sort(this->send_list.begin(), 
00715             this->send_list.end());
00716   
00717   // Now use std::unique to remove duplicate entries
00718   std::vector<unsigned int>::iterator new_end =
00719     std::unique (this->send_list.begin(), 
00720                  this->send_list.end());
00721   
00722   // Remove the end of the send_list.  Use the "swap trick"
00723   // from Effective STL
00724   std::vector<unsigned int> 
00725     (this->send_list.begin(), new_end).swap (this->send_list);
00726 }
00727 
00728 
00729 
00730 #ifndef LIBMESH_ENABLE_AMR
00731 void System::BuildProjectionList::operator()(const ConstElemRange &)
00732 {
00733   libmesh_error();
00734 #else
00735 void System::BuildProjectionList::operator()(const ConstElemRange &range)
00736 {
00737   // The DofMap for this system
00738   const DofMap& dof_map = system.get_dof_map();
00739 
00740   // We can handle all the variables at once.
00741   // The old global DOF indices
00742   std::vector<unsigned int> di;
00743    
00744   // Iterate over the elements in the range
00745   for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it)
00746     {
00747       const Elem* elem = *elem_it;
00748       const Elem* parent = elem->parent();
00749       
00750       if (elem->refinement_flag() == Elem::JUST_REFINED)
00751         {
00752           libmesh_assert (parent != NULL);
00753           unsigned int old_parent_level = parent->p_level();
00754           
00755           if (elem->p_refinement_flag() == Elem::JUST_REFINED)
00756             {
00757               // We may have done p refinement or coarsening as well;
00758               // if so then we need to reset the parent's p level
00759               // so we can get the right DoFs from it
00760               libmesh_assert(elem->p_level() > 0);
00761               (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() - 1);
00762             }
00763           else if (elem->p_refinement_flag() == Elem::JUST_COARSENED)
00764             {
00765               (const_cast<Elem *>(parent))->hack_p_level(elem->p_level() + 1);
00766             }
00767 
00768           dof_map.old_dof_indices (parent, di);
00769           
00770           // Fix up the parent's p level in case we changed it
00771           (const_cast<Elem *>(parent))->hack_p_level(old_parent_level);
00772         }
00773       else
00774         dof_map.old_dof_indices (elem, di);
00775       
00776       this->send_list.insert(send_list.end(), di.begin(), di.end());
00777     }  // end elem loop
00778 #endif // LIBMESH_ENABLE_AMR
00779 }
00780 
00781 
00782 
00783 void System::BuildProjectionList::join(const BuildProjectionList &other)
00784 {
00785   // Joining simply requires I add the dof indices from the other object
00786   this->send_list.insert(this->send_list.end(),
00787                          other.send_list.begin(),
00788                          other.send_list.end());
00789 }
00790 
00791 
00792 void System::ProjectSolution::operator()(const ConstElemRange &range) const
00793 {
00794   // We need data to project
00795   libmesh_assert(fptr);
00796 
00804   // The number of variables in this system
00805   const unsigned int n_variables = system.n_vars();
00806 
00807   // The dimensionality of the current mesh
00808   const unsigned int dim = system.get_mesh().mesh_dimension();
00809 
00810   // The DofMap for this system
00811   const DofMap& dof_map = system.get_dof_map();
00812 
00813   // The element matrix and RHS for projections.
00814   // Note that Ke is always real-valued, whereas
00815   // Fe may be complex valued if complex number
00816   // support is enabled
00817   DenseMatrix<Real> Ke;
00818   DenseVector<Number> Fe;
00819   // The new element coefficients
00820   DenseVector<Number> Ue;
00821 
00822 
00823   // Loop over all the variables in the system
00824   for (unsigned int var=0; var<n_variables; var++)
00825     {
00826       if (dof_map.variable(var).type().family == SCALAR)
00827         continue;
00828 
00829       // Get FE objects of the appropriate type
00830       const FEType& fe_type = dof_map.variable_type(var);     
00831       AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));      
00832 
00833       // Prepare variables for projection
00834       AutoPtr<QBase> qrule     (fe_type.default_quadrature_rule(dim));
00835       AutoPtr<QBase> qedgerule (fe_type.default_quadrature_rule(1));
00836       AutoPtr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1));
00837 
00838       // The values of the shape functions at the quadrature
00839       // points
00840       const std::vector<std::vector<Real> >& phi = fe->get_phi();
00841 
00842       // The gradients of the shape functions at the quadrature
00843       // points on the child element.
00844       const std::vector<std::vector<RealGradient> > *dphi = NULL;
00845 
00846       const FEContinuity cont = fe->get_continuity();
00847 
00848       if (cont == C_ONE)
00849         {
00850           // We'll need gradient data for a C1 projection
00851           libmesh_assert(gptr);
00852 
00853           const std::vector<std::vector<RealGradient> >&
00854             ref_dphi = fe->get_dphi();
00855           dphi = &ref_dphi;
00856         }
00857 
00858       // The Jacobian * quadrature weight at the quadrature points
00859       const std::vector<Real>& JxW =
00860         fe->get_JxW();
00861      
00862       // The XYZ locations of the quadrature points
00863       const std::vector<Point>& xyz_values =
00864         fe->get_xyz();
00865 
00866       // The global DOF indices
00867       std::vector<unsigned int> dof_indices;
00868       // Side/edge DOF indices
00869       std::vector<unsigned int> side_dofs;
00870    
00871       // Iterate over all the elements in the range
00872       for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it)
00873         {
00874           const Elem* elem = *elem_it;
00875 
00876           // Update the DOF indices for this element based on
00877           // the current mesh
00878           dof_map.dof_indices (elem, dof_indices, var);
00879 
00880           // The number of DOFs on the element
00881           const unsigned int n_dofs = dof_indices.size();
00882 
00883           // Fixed vs. free DoFs on edge/face projections
00884           std::vector<char> dof_is_fixed(n_dofs, false); // bools
00885           std::vector<int> free_dof(n_dofs, 0);
00886 
00887           // The element type
00888           const ElemType elem_type = elem->type();
00889 
00890           // The number of nodes on the new element
00891           const unsigned int n_nodes = elem->n_nodes();
00892 
00893           // Zero the interpolated values
00894           Ue.resize (n_dofs); Ue.zero();
00895 
00896           // In general, we need a series of
00897           // projections to ensure a unique and continuous
00898           // solution.  We start by interpolating nodes, then
00899           // hold those fixed and project edges, then
00900           // hold those fixed and project faces, then
00901           // hold those fixed and project interiors
00902 
00903           // Interpolate node values first
00904           unsigned int current_dof = 0;
00905           for (unsigned int n=0; n!= n_nodes; ++n)
00906             {
00907               // FIXME: this should go through the DofMap,
00908               // not duplicate dof_indices code badly!
00909               const unsigned int nc =
00910                 FEInterface::n_dofs_at_node (dim, fe_type, elem_type,
00911                                              n);
00912               if (!elem->is_vertex(n))
00913                 {
00914                   current_dof += nc;
00915                   continue;
00916                 }
00917               if (cont == DISCONTINUOUS)
00918                 {
00919                   libmesh_assert(nc == 0);
00920                 }
00921               // Assume that C_ZERO elements have a single nodal
00922               // value shape function
00923               else if (cont == C_ZERO)
00924                 {
00925                   libmesh_assert(nc == 1);
00926                   Ue(current_dof) = fptr(elem->point(n),
00927                                          parameters,
00928                                          system.name(),
00929                                          system.variable_name(var));
00930                   dof_is_fixed[current_dof] = true;
00931                   current_dof++;
00932                 }
00933               // The hermite element vertex shape functions are weird
00934               else if (fe_type.family == HERMITE)
00935                 {
00936                   Ue(current_dof) = fptr(elem->point(n),
00937                                          parameters,
00938                                          system.name(),
00939                                          system.variable_name(var));
00940                   dof_is_fixed[current_dof] = true;
00941                   current_dof++;
00942                   Gradient g = gptr(elem->point(n),
00943                                     parameters,
00944                                     system.name(),
00945                                     system.variable_name(var));
00946                   // x derivative
00947                   Ue(current_dof) = g(0);
00948                   dof_is_fixed[current_dof] = true;
00949                   current_dof++;
00950                   if (dim > 1)
00951                     {
00952                       // We'll finite difference mixed derivatives
00953                       Point nxminus = elem->point(n),
00954                             nxplus = elem->point(n);
00955                       nxminus(0) -= TOLERANCE;
00956                       nxplus(0) += TOLERANCE;
00957                       Gradient gxminus = gptr(nxminus,
00958                                               parameters,
00959                                               system.name(),
00960                                               system.variable_name(var));
00961                       Gradient gxplus = gptr(nxminus,
00962                                              parameters,
00963                                              system.name(),
00964                                              system.variable_name(var));
00965                       // y derivative
00966                       Ue(current_dof) = g(1);
00967                       dof_is_fixed[current_dof] = true;
00968                       current_dof++;
00969                       // xy derivative
00970                       Ue(current_dof) = (gxplus(1) - gxminus(1))
00971                                         / 2. / TOLERANCE;
00972                       dof_is_fixed[current_dof] = true;
00973                       current_dof++;
00974 
00975                       if (dim > 2)
00976                         {
00977                           // z derivative
00978                           Ue(current_dof) = g(2);
00979                           dof_is_fixed[current_dof] = true;
00980                           current_dof++;
00981                           // xz derivative
00982                           Ue(current_dof) = (gxplus(2) - gxminus(2))
00983                                             / 2. / TOLERANCE;
00984                           dof_is_fixed[current_dof] = true;
00985                           current_dof++;
00986                           // We need new points for yz
00987                           Point nyminus = elem->point(n),
00988                                 nyplus = elem->point(n);
00989                           nyminus(1) -= TOLERANCE;
00990                           nyplus(1) += TOLERANCE;
00991                           Gradient gyminus = gptr(nyminus,
00992                                                   parameters,
00993                                                   system.name(),
00994                                                   system.variable_name(var));
00995                           Gradient gyplus = gptr(nyminus,
00996                                                  parameters,
00997                                                  system.name(),
00998                                                  system.variable_name(var));
00999                           // xz derivative
01000                           Ue(current_dof) = (gyplus(2) - gyminus(2))
01001                                             / 2. / TOLERANCE;
01002                           dof_is_fixed[current_dof] = true;
01003                           current_dof++;
01004                           // Getting a 2nd order xyz is more tedious
01005                           Point nxmym = elem->point(n),
01006                                 nxmyp = elem->point(n),
01007                                 nxpym = elem->point(n),
01008                                 nxpyp = elem->point(n);
01009                           nxmym(0) -= TOLERANCE;
01010                           nxmym(1) -= TOLERANCE;
01011                           nxmyp(0) -= TOLERANCE;
01012                           nxmyp(1) += TOLERANCE;
01013                           nxpym(0) += TOLERANCE;
01014                           nxpym(1) -= TOLERANCE;
01015                           nxpyp(0) += TOLERANCE;
01016                           nxpyp(1) += TOLERANCE;
01017                           Gradient gxmym = gptr(nxmym,
01018                                                 parameters,
01019                                                 system.name(),
01020                                                 system.variable_name(var));
01021                           Gradient gxmyp = gptr(nxmyp,
01022                                                 parameters,
01023                                                 system.name(),
01024                                                 system.variable_name(var));
01025                           Gradient gxpym = gptr(nxpym,
01026                                                 parameters,
01027                                                 system.name(),
01028                                                 system.variable_name(var));
01029                           Gradient gxpyp = gptr(nxpyp,
01030                                                 parameters,
01031                                                 system.name(),
01032                                                 system.variable_name(var));
01033                           Number gxzplus = (gxpyp(2) - gxmyp(2))
01034                                          / 2. / TOLERANCE;
01035                           Number gxzminus = (gxpym(2) - gxmym(2))
01036                                           / 2. / TOLERANCE;
01037                           // xyz derivative
01038                           Ue(current_dof) = (gxzplus - gxzminus)
01039                                             / 2. / TOLERANCE;
01040                           dof_is_fixed[current_dof] = true;
01041                           current_dof++;
01042                         }
01043                     }
01044                 }
01045               // Assume that other C_ONE elements have a single nodal
01046               // value shape function and nodal gradient component
01047               // shape functions
01048               else if (cont == C_ONE)
01049                 {
01050                   libmesh_assert(nc == 1 + dim);
01051                   Ue(current_dof) = fptr(elem->point(n),
01052                                          parameters,
01053                                          system.name(),
01054                                          system.variable_name(var));
01055                   dof_is_fixed[current_dof] = true;
01056                   current_dof++;
01057                   Gradient g = gptr(elem->point(n),
01058                                     parameters,
01059                                     system.name(),
01060                                     system.variable_name(var));
01061                   for (unsigned int i=0; i!= dim; ++i)
01062                     {
01063                       Ue(current_dof) = g(i);
01064                       dof_is_fixed[current_dof] = true;
01065                       current_dof++;
01066                     }
01067                 }
01068               else
01069                 libmesh_error();
01070             }
01071 
01072           // In 3D, project any edge values next
01073           if (dim > 2 && cont != DISCONTINUOUS)
01074             for (unsigned int e=0; e != elem->n_edges(); ++e)
01075               {
01076                 FEInterface::dofs_on_edge(elem, dim, fe_type, e,
01077                                           side_dofs);
01078 
01079                 // Some edge dofs are on nodes and already
01080                 // fixed, others are free to calculate
01081                 unsigned int free_dofs = 0;
01082                 for (unsigned int i=0; i != side_dofs.size(); ++i)
01083                   if (!dof_is_fixed[side_dofs[i]])
01084                     free_dof[free_dofs++] = i;
01085 
01086                 // There may be nothing to project
01087                 if (!free_dofs)
01088                   continue;
01089 
01090                 Ke.resize (free_dofs, free_dofs); Ke.zero();
01091                 Fe.resize (free_dofs); Fe.zero();
01092                 // The new edge coefficients
01093                 DenseVector<Number> Uedge(free_dofs);
01094 
01095                 // Initialize FE data on the edge
01096                 fe->attach_quadrature_rule (qedgerule.get());
01097                 fe->edge_reinit (elem, e);
01098                 const unsigned int n_qp = qedgerule->n_points();
01099 
01100                 // Loop over the quadrature points
01101                 for (unsigned int qp=0; qp<n_qp; qp++)
01102                   {
01103                     // solution at the quadrature point       
01104                     Number fineval = fptr(xyz_values[qp],
01105                                           parameters,
01106                                           system.name(),
01107                                           system.variable_name(var));
01108                     // solution grad at the quadrature point          
01109                     Gradient finegrad;
01110                     if (cont == C_ONE)
01111                       finegrad = gptr(xyz_values[qp], parameters,
01112                                       system.name(),
01113                                       system.variable_name(var));
01114 
01115                     // Form edge projection matrix
01116                     for (unsigned int sidei=0, freei=0; 
01117                          sidei != side_dofs.size(); ++sidei)
01118                       {
01119                         unsigned int i = side_dofs[sidei];
01120                         // fixed DoFs aren't test functions
01121                         if (dof_is_fixed[i])
01122                           continue;
01123                         for (unsigned int sidej=0, freej=0;
01124                              sidej != side_dofs.size(); ++sidej)
01125                           {
01126                             unsigned int j = side_dofs[sidej];
01127                             if (dof_is_fixed[j])
01128                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
01129                                            JxW[qp] * Ue(j);
01130                             else
01131                               Ke(freei,freej) += phi[i][qp] *
01132                                                  phi[j][qp] * JxW[qp];
01133                             if (cont == C_ONE)
01134                               {
01135                                 if (dof_is_fixed[j])
01136                                   Fe(freei) -= ((*dphi)[i][qp] *
01137                                                 (*dphi)[j][qp]) *
01138                                                 JxW[qp] * Ue(j);
01139                                 else
01140                                   Ke(freei,freej) += ((*dphi)[i][qp] *
01141                                                       (*dphi)[j][qp])
01142                                                       * JxW[qp];
01143                               }
01144                             if (!dof_is_fixed[j])
01145                               freej++;
01146                           }
01147                         Fe(freei) += phi[i][qp] * fineval * JxW[qp];
01148                         if (cont == C_ONE)
01149                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
01150                                        JxW[qp];
01151                         freei++;
01152                       }
01153                   }
01154 
01155                 Ke.cholesky_solve(Fe, Uedge);
01156 
01157                 // Transfer new edge solutions to element
01158                 for (unsigned int i=0; i != free_dofs; ++i)
01159                   {
01160                     Number &ui = Ue(side_dofs[free_dof[i]]);
01161                     libmesh_assert(std::abs(ui) < TOLERANCE ||
01162                            std::abs(ui - Uedge(i)) < TOLERANCE);
01163                     ui = Uedge(i);
01164                     dof_is_fixed[side_dofs[free_dof[i]]] = true;
01165                   }
01166               }
01167                  
01168           // Project any side values (edges in 2D, faces in 3D)
01169           if (dim > 1 && cont != DISCONTINUOUS)
01170             for (unsigned int s=0; s != elem->n_sides(); ++s)
01171               {
01172                 FEInterface::dofs_on_side(elem, dim, fe_type, s,
01173                                           side_dofs);
01174 
01175                 // Some side dofs are on nodes/edges and already
01176                 // fixed, others are free to calculate
01177                 unsigned int free_dofs = 0;
01178                 for (unsigned int i=0; i != side_dofs.size(); ++i)
01179                   if (!dof_is_fixed[side_dofs[i]])
01180                     free_dof[free_dofs++] = i;
01181 
01182                 // There may be nothing to project
01183                 if (!free_dofs)
01184                   continue;
01185              
01186                 Ke.resize (free_dofs, free_dofs); Ke.zero();
01187                 Fe.resize (free_dofs); Fe.zero();
01188                 // The new side coefficients
01189                 DenseVector<Number> Uside(free_dofs);
01190 
01191                 // Initialize FE data on the side
01192                 fe->attach_quadrature_rule (qsiderule.get());
01193                 fe->reinit (elem, s);
01194                 const unsigned int n_qp = qsiderule->n_points();
01195 
01196                 // Loop over the quadrature points
01197                 for (unsigned int qp=0; qp<n_qp; qp++)
01198                   {
01199                     // solution at the quadrature point       
01200                     Number fineval = fptr(xyz_values[qp],
01201                                           parameters,
01202                                           system.name(),
01203                                           system.variable_name(var));
01204                     // solution grad at the quadrature point          
01205                     Gradient finegrad;
01206                     if (cont == C_ONE)
01207                       finegrad = gptr(xyz_values[qp], parameters,
01208                                       system.name(),
01209                                       system.variable_name(var));
01210 
01211                     // Form side projection matrix
01212                     for (unsigned int sidei=0, freei=0;
01213                          sidei != side_dofs.size(); ++sidei)
01214                       {
01215                         unsigned int i = side_dofs[sidei];
01216                         // fixed DoFs aren't test functions
01217                         if (dof_is_fixed[i])
01218                           continue;
01219                         for (unsigned int sidej=0, freej=0;
01220                              sidej != side_dofs.size(); ++sidej)
01221                           {
01222                             unsigned int j = side_dofs[sidej];
01223                             if (dof_is_fixed[j])
01224                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
01225                                            JxW[qp] * Ue(j);
01226                             else
01227                               Ke(freei,freej) += phi[i][qp] *
01228                                                  phi[j][qp] * JxW[qp];
01229                             if (cont == C_ONE)
01230                               {
01231                                 if (dof_is_fixed[j])
01232                                   Fe(freei) -= ((*dphi)[i][qp] *
01233                                                 (*dphi)[j][qp]) *
01234                                                JxW[qp] * Ue(j);
01235                                 else
01236                                   Ke(freei,freej) += ((*dphi)[i][qp] *
01237                                                       (*dphi)[j][qp])
01238                                                      * JxW[qp];
01239                               }
01240                             if (!dof_is_fixed[j])
01241                               freej++;
01242                           }
01243                         Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
01244                         if (cont == C_ONE)
01245                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
01246                                        JxW[qp];
01247                         freei++;
01248                       }
01249                   }
01250 
01251                 Ke.cholesky_solve(Fe, Uside);
01252 
01253                 // Transfer new side solutions to element
01254                 for (unsigned int i=0; i != free_dofs; ++i)
01255                   {
01256                     Number &ui = Ue(side_dofs[free_dof[i]]);
01257                     libmesh_assert(std::abs(ui) < TOLERANCE ||
01258                            std::abs(ui - Uside(i)) < TOLERANCE);
01259                     ui = Uside(i);
01260                     dof_is_fixed[side_dofs[free_dof[i]]] = true;
01261                   }
01262               }
01263 
01264           // Project the interior values, finally
01265 
01266           // Some interior dofs are on nodes/edges/sides and
01267           // already fixed, others are free to calculate
01268           unsigned int free_dofs = 0;
01269           for (unsigned int i=0; i != n_dofs; ++i)
01270             if (!dof_is_fixed[i])
01271               free_dof[free_dofs++] = i;
01272 
01273           // There may be nothing to project
01274           if (free_dofs)
01275             {
01276              
01277           Ke.resize (free_dofs, free_dofs); Ke.zero();
01278           Fe.resize (free_dofs); Fe.zero();
01279           // The new interior coefficients
01280           DenseVector<Number> Uint(free_dofs);
01281 
01282           // Initialize FE data
01283           fe->attach_quadrature_rule (qrule.get());
01284           fe->reinit (elem);
01285           const unsigned int n_qp = qrule->n_points();
01286 
01287           // Loop over the quadrature points
01288           for (unsigned int qp=0; qp<n_qp; qp++)
01289             {
01290               // solution at the quadrature point             
01291               Number fineval = fptr(xyz_values[qp],
01292                                     parameters,
01293                                     system.name(),
01294                                     system.variable_name(var));
01295               // solution grad at the quadrature point        
01296               Gradient finegrad;
01297               if (cont == C_ONE)
01298                 finegrad = gptr(xyz_values[qp], parameters,
01299                                 system.name(),
01300                                 system.variable_name(var));
01301 
01302               // Form interior projection matrix
01303               for (unsigned int i=0, freei=0; i != n_dofs; ++i)
01304                 {
01305                   // fixed DoFs aren't test functions
01306                   if (dof_is_fixed[i])
01307                     continue;
01308                   for (unsigned int j=0, freej=0; j != n_dofs; ++j)
01309                     {
01310                       if (dof_is_fixed[j])
01311                         Fe(freei) -= phi[i][qp] * phi[j][qp] * JxW[qp]
01312                                      * Ue(j);
01313                       else
01314                         Ke(freei,freej) += phi[i][qp] * phi[j][qp] *
01315                                            JxW[qp];
01316                       if (cont == C_ONE)
01317                         {
01318                           if (dof_is_fixed[j])
01319                             Fe(freei) -= ((*dphi)[i][qp] *
01320                                          (*dphi)[j][qp]) * JxW[qp] *
01321                                          Ue(j);
01322                           else
01323                             Ke(freei,freej) += ((*dphi)[i][qp] *
01324                                                 (*dphi)[j][qp]) *
01325                                                JxW[qp];
01326                         }
01327                       if (!dof_is_fixed[j])
01328                         freej++;
01329                     }
01330                   Fe(freei) += phi[i][qp] * fineval * JxW[qp];
01331                   if (cont == C_ONE)
01332                     Fe(freei) += (finegrad * (*dphi)[i][qp]) * JxW[qp];
01333                   freei++;
01334                 }
01335             }
01336           Ke.cholesky_solve(Fe, Uint);
01337 
01338           // Transfer new interior solutions to element
01339           for (unsigned int i=0; i != free_dofs; ++i)
01340             {
01341               Number &ui = Ue(free_dof[i]);
01342               libmesh_assert(std::abs(ui) < TOLERANCE ||
01343                      std::abs(ui - Uint(i)) < TOLERANCE);
01344               ui = Uint(i);
01345               dof_is_fixed[free_dof[i]] = true;
01346             }
01347 
01348             } // if there are free interior dofs
01349 
01350           // Make sure every DoF got reached!
01351           for (unsigned int i=0; i != n_dofs; ++i)
01352             libmesh_assert(dof_is_fixed[i]);
01353 
01354           const unsigned int
01355             first = new_vector.first_local_index(),
01356             last  = new_vector.last_local_index();
01357           
01358           // Lock the new_vector since it is shared among threads.
01359           {
01360             Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
01361             
01362             for (unsigned int i = 0; i < n_dofs; i++) 
01363               // We may be projecting a new zero value onto
01364               // an old nonzero approximation - RHS
01365               // if (Ue(i) != 0.)
01366               if ((dof_indices[i] >= first) &&
01367                   (dof_indices[i] <  last))
01368                 {
01369                   new_vector.set(dof_indices[i], Ue(i));
01370                 }
01371           }
01372         }  // end elem loop
01373     } // end variables loop
01374 }

Site Created By: libMesh Developers
Last modified: November 25 2009 03:43:48.

Hosted By:
SourceForge.net Logo