libMesh::ProjectSolution Class Reference

List of all members.

Public Member Functions

 ProjectSolution (const System &system_in, FunctionBase< Number > *f_in, FunctionBase< Gradient > *g_in, const Parameters &parameters_in, NumericVector< Number > &new_v_in)
 ProjectSolution (const ProjectSolution &in)
void operator() (const ConstElemRange &range) const

Private Attributes

const Systemsystem
AutoPtr< FunctionBase< Number > > f
AutoPtr< FunctionBase< Gradient > > g
const Parametersparameters
NumericVector< Number > & new_vector

Detailed Description

This class implements projecting an arbitrary function to the current mesh. This may be exectued in parallel on multiple threads.

Definition at line 110 of file system_projection.C.


Constructor & Destructor Documentation

libMesh::ProjectSolution::ProjectSolution ( const System system_in,
FunctionBase< Number > *  f_in,
FunctionBase< Gradient > *  g_in,
const Parameters parameters_in,
NumericVector< Number > &  new_v_in 
) [inline]

Definition at line 121 of file system_projection.C.

References f, g, and libMesh::AutoPtr< Tp >::get().

00125                                                       :
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     }

libMesh::ProjectSolution::ProjectSolution ( const ProjectSolution in  )  [inline]

Definition at line 138 of file system_projection.C.

References f, g, and libMesh::AutoPtr< Tp >::get().

00138                                                 :
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     }


Member Function Documentation

void libMesh::ProjectSolution::operator() ( const ConstElemRange range  )  const

This method projects an arbitrary solution to the current mesh. The input function f gives the arbitrary solution, while the new_vector (which should already be correctly sized) gives the solution (to be computed) on the current mesh.

Definition at line 1250 of file system_projection.C.

References std::abs(), libMesh::Variable::active_on_subdomain(), libMesh::StoredRange< iterator_type, object_type >::begin(), libMesh::FEGenericBase< Real >::build(), libMeshEnums::C_ONE, libMeshEnums::C_ZERO, libMesh::DenseMatrix< T >::cholesky_solve(), libMesh::FunctionBase< Output >::component(), libMesh::FEType::default_quadrature_rule(), libMeshEnums::DISCONTINUOUS, libMesh::DofMap::dof_indices(), libMesh::FEInterface::dofs_on_edge(), libMesh::FEInterface::dofs_on_side(), libMesh::StoredRange< iterator_type, object_type >::end(), libMesh::FEType::family, libMesh::NumericVector< T >::first_local_index(), libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMeshEnums::HERMITE, libMesh::Elem::is_vertex(), libMesh::NumericVector< T >::last_local_index(), libMesh::MeshBase::mesh_dimension(), libMesh::FEInterface::n_dofs_at_node(), libMesh::Elem::n_edges(), libMesh::Elem::n_nodes(), libMesh::Elem::n_sides(), libMesh::System::n_vars(), libMesh::Elem::point(), libMesh::DenseMatrix< T >::resize(), libMesh::DenseVector< T >::resize(), libMeshEnums::SCALAR, libMesh::NumericVector< T >::set(), libMesh::Threads::spin_mtx, libMesh::Elem::subdomain_id(), system, libMesh::System::time, libMesh::TOLERANCE, libMesh::Elem::type(), libMesh::Variable::type(), libMesh::DofMap::variable(), libMesh::System::variable_scalar_number(), libMesh::DenseMatrix< T >::zero(), and libMesh::DenseVector< T >::zero().

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 }


Member Data Documentation

Definition at line 115 of file system_projection.C.

Referenced by ProjectSolution().

Definition at line 116 of file system_projection.C.

Referenced by ProjectSolution().

Definition at line 117 of file system_projection.C.

Definition at line 113 of file system_projection.C.

Referenced by operator()().


The documentation for this class was generated from the following file:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:06 UTC

Hosted By:
SourceForge.net Logo