libMesh::BoundaryProjectSolution Class Reference

List of all members.

Public Member Functions

 BoundaryProjectSolution (const std::set< boundary_id_type > &b_in, const std::vector< unsigned int > &variables_in, const System &system_in, FunctionBase< Number > *f_in, FunctionBase< Gradient > *g_in, const Parameters &parameters_in, NumericVector< Number > &new_v_in)
 BoundaryProjectSolution (const BoundaryProjectSolution &in)
void operator() (const ConstElemRange &range) const

Private Attributes

const std::set
< boundary_id_type > & 
b
const std::vector< unsigned int > & variables
const Systemsystem
AutoPtr< FunctionBase< Number > > f
AutoPtr< FunctionBase< Gradient > > g
const Parametersparameters
NumericVector< Number > & new_vector

Detailed Description

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

Definition at line 200 of file system_projection.C.


Constructor & Destructor Documentation

libMesh::BoundaryProjectSolution::BoundaryProjectSolution ( const std::set< boundary_id_type > &  b_in,
const std::vector< unsigned int > &  variables_in,
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 212 of file system_projection.C.

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

00218                                                               :
00219     b(b_in),
00220     variables(variables_in),
00221     system(system_in),
00222     f(f_in ? f_in->clone() : AutoPtr<FunctionBase<Number> >(NULL)),
00223     g(g_in ? g_in->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)),
00224     parameters(parameters_in),
00225     new_vector(new_v_in)
00226     {
00227       libmesh_assert(f.get());
00228       f->init();
00229       if (g.get())
00230         g->init();
00231     }

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

Definition at line 233 of file system_projection.C.

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

00233                                                                 :
00234     b(in.b),
00235     variables(in.variables),
00236     system(in.system),
00237     f(in.f.get() ? in.f->clone() : AutoPtr<FunctionBase<Number> >(NULL)),
00238     g(in.g.get() ? in.g->clone() : AutoPtr<FunctionBase<Gradient> >(NULL)),
00239     parameters(in.parameters),
00240     new_vector(in.new_vector)
00241     {
00242       libmesh_assert(f.get());
00243       f->init();
00244       if (g.get())
00245         g->init();
00246     }


Member Function Documentation

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

This method projects an arbitrary boundary 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 2463 of file system_projection.C.

References std::abs(), libMesh::Variable::active_on_subdomain(), libMesh::StoredRange< iterator_type, object_type >::begin(), libMesh::BoundaryInfo::boundary_ids(), libMesh::MeshBase::boundary_info, 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_edge_on_side(), libMesh::Elem::is_node_on_side(), 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::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().

02464 {
02465   // We need data to project
02466   libmesh_assert(f.get());
02467 
02475   // The dimensionality of the current mesh
02476   const unsigned int dim = system.get_mesh().mesh_dimension();
02477 
02478   // The DofMap for this system
02479   const DofMap& dof_map = system.get_dof_map();
02480 
02481   // Boundary info for the current mesh
02482   const BoundaryInfo& boundary_info = *system.get_mesh().boundary_info;
02483 
02484   // The element matrix and RHS for projections.
02485   // Note that Ke is always real-valued, whereas
02486   // Fe may be complex valued if complex number
02487   // support is enabled
02488   DenseMatrix<Real> Ke;
02489   DenseVector<Number> Fe;
02490   // The new element coefficients
02491   DenseVector<Number> Ue;
02492 
02493 
02494   // Loop over all the variables we've been requested to project
02495   for (unsigned int v=0; v!=variables.size(); v++)
02496     {
02497       const unsigned int var = variables[v];
02498 
02499       const Variable& variable = dof_map.variable(var);
02500 
02501       const FEType& fe_type = variable.type();
02502 
02503       if (fe_type.family == SCALAR)
02504         continue;
02505 
02506       const unsigned int var_component =
02507         system.variable_scalar_number(var, 0);
02508 
02509       // Get FE objects of the appropriate type
02510       AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));
02511 
02512       // Prepare variables for projection
02513       AutoPtr<QBase> qedgerule (fe_type.default_quadrature_rule(1));
02514       AutoPtr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1));
02515 
02516       // The values of the shape functions at the quadrature
02517       // points
02518       const std::vector<std::vector<Real> >& phi = fe->get_phi();
02519 
02520       // The gradients of the shape functions at the quadrature
02521       // points on the child element.
02522       const std::vector<std::vector<RealGradient> > *dphi = NULL;
02523 
02524       const FEContinuity cont = fe->get_continuity();
02525 
02526       if (cont == C_ONE)
02527         {
02528           // We'll need gradient data for a C1 projection
02529           libmesh_assert(g.get());
02530 
02531           const std::vector<std::vector<RealGradient> >&
02532             ref_dphi = fe->get_dphi();
02533           dphi = &ref_dphi;
02534         }
02535 
02536       // The Jacobian * quadrature weight at the quadrature points
02537       const std::vector<Real>& JxW =
02538         fe->get_JxW();
02539 
02540       // The XYZ locations of the quadrature points
02541       const std::vector<Point>& xyz_values =
02542         fe->get_xyz();
02543 
02544       // The global DOF indices
02545       std::vector<dof_id_type> dof_indices;
02546       // Side/edge DOF indices
02547       std::vector<unsigned int> side_dofs;
02548 
02549       // Iterate over all the elements in the range
02550       for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it != range.end(); ++elem_it)
02551         {
02552           const Elem* elem = *elem_it;
02553 
02554           // Per-subdomain variables don't need to be projected on
02555           // elements where they're not active
02556           if (!variable.active_on_subdomain(elem->subdomain_id()))
02557             continue;
02558 
02559           // Find out which nodes, edges and sides are on a requested
02560           // boundary:
02561           std::vector<bool> is_boundary_node(elem->n_nodes(), false),
02562                             is_boundary_edge(elem->n_edges(), false),
02563                             is_boundary_side(elem->n_sides(), false);
02564           for (unsigned char s=0; s != elem->n_sides(); ++s)
02565             {
02566               // First see if this side has been requested
02567               const std::vector<boundary_id_type>& bc_ids =
02568                 boundary_info.boundary_ids (elem, s);
02569               bool do_this_side = false;
02570               for (unsigned int i=0; i != bc_ids.size(); ++i)
02571                 if (b.count(bc_ids[i]))
02572                   {
02573                     do_this_side = true;
02574                     break;
02575                   }
02576               if (!do_this_side)
02577                 continue;
02578 
02579               is_boundary_side[s] = true;
02580 
02581               // Then see what nodes and what edges are on it
02582               for (unsigned int n=0; n != elem->n_nodes(); ++n)
02583                 if (elem->is_node_on_side(n,s))
02584                   is_boundary_node[n] = true;
02585               for (unsigned int e=0; e != elem->n_edges(); ++e)
02586                 if (elem->is_edge_on_side(e,s))
02587                   is_boundary_edge[e] = true;
02588             }
02589 
02590           // Update the DOF indices for this element based on
02591           // the current mesh
02592           dof_map.dof_indices (elem, dof_indices, var);
02593 
02594           // The number of DOFs on the element
02595           const unsigned int n_dofs =
02596             libmesh_cast_int<unsigned int>(dof_indices.size());
02597 
02598           // Fixed vs. free DoFs on edge/face projections
02599           std::vector<char> dof_is_fixed(n_dofs, false); // bools
02600           std::vector<int> free_dof(n_dofs, 0);
02601 
02602           // The element type
02603           const ElemType elem_type = elem->type();
02604 
02605           // The number of nodes on the new element
02606           const unsigned int n_nodes = elem->n_nodes();
02607 
02608           // Zero the interpolated values
02609           Ue.resize (n_dofs); Ue.zero();
02610 
02611           // In general, we need a series of
02612           // projections to ensure a unique and continuous
02613           // solution.  We start by interpolating boundary nodes, then
02614           // hold those fixed and project boundary edges, then hold
02615           // those fixed and project boundary faces,
02616 
02617           // Interpolate node values first
02618           unsigned int current_dof = 0;
02619           for (unsigned int n=0; n!= n_nodes; ++n)
02620             {
02621               // FIXME: this should go through the DofMap,
02622               // not duplicate dof_indices code badly!
02623               const unsigned int nc =
02624                 FEInterface::n_dofs_at_node (dim, fe_type, elem_type,
02625                                              n);
02626               if (!elem->is_vertex(n) || !is_boundary_node[n])
02627                 {
02628                   current_dof += nc;
02629                   continue;
02630                 }
02631               if (cont == DISCONTINUOUS)
02632                 {
02633                   libmesh_assert_equal_to (nc, 0);
02634                 }
02635               // Assume that C_ZERO elements have a single nodal
02636               // value shape function
02637               else if (cont == C_ZERO)
02638                 {
02639                   libmesh_assert_equal_to (nc, 1);
02640                   Ue(current_dof) = f->component(var_component,
02641                                                  elem->point(n),
02642                                                  system.time);
02643                   dof_is_fixed[current_dof] = true;
02644                   current_dof++;
02645                 }
02646               // The hermite element vertex shape functions are weird
02647               else if (fe_type.family == HERMITE)
02648                 {
02649                   Ue(current_dof) = f->component(var_component,
02650                                                  elem->point(n),
02651                                                  system.time);
02652                   dof_is_fixed[current_dof] = true;
02653                   current_dof++;
02654                   Gradient grad = g->component(var_component,
02655                                                elem->point(n),
02656                                                system.time);
02657                   // x derivative
02658                   Ue(current_dof) = grad(0);
02659                   dof_is_fixed[current_dof] = true;
02660                   current_dof++;
02661                   if (dim > 1)
02662                     {
02663                       // We'll finite difference mixed derivatives
02664                       Point nxminus = elem->point(n),
02665                             nxplus = elem->point(n);
02666                       nxminus(0) -= TOLERANCE;
02667                       nxplus(0) += TOLERANCE;
02668                       Gradient gxminus = g->component(var_component,
02669                                                       nxminus,
02670                                                       system.time);
02671                       Gradient gxplus = g->component(var_component,
02672                                                      nxplus,
02673                                                      system.time);
02674                       // y derivative
02675                       Ue(current_dof) = grad(1);
02676                       dof_is_fixed[current_dof] = true;
02677                       current_dof++;
02678                       // xy derivative
02679                       Ue(current_dof) = (gxplus(1) - gxminus(1))
02680                                         / 2. / TOLERANCE;
02681                       dof_is_fixed[current_dof] = true;
02682                       current_dof++;
02683 
02684                       if (dim > 2)
02685                         {
02686                           // z derivative
02687                           Ue(current_dof) = grad(2);
02688                           dof_is_fixed[current_dof] = true;
02689                           current_dof++;
02690                           // xz derivative
02691                           Ue(current_dof) = (gxplus(2) - gxminus(2))
02692                                             / 2. / TOLERANCE;
02693                           dof_is_fixed[current_dof] = true;
02694                           current_dof++;
02695                           // We need new points for yz
02696                           Point nyminus = elem->point(n),
02697                                 nyplus = elem->point(n);
02698                           nyminus(1) -= TOLERANCE;
02699                           nyplus(1) += TOLERANCE;
02700                           Gradient gyminus = g->component(var_component,
02701                                                           nyminus,
02702                                                           system.time);
02703                           Gradient gyplus = g->component(var_component,
02704                                                          nyplus,
02705                                                          system.time);
02706                           // xz derivative
02707                           Ue(current_dof) = (gyplus(2) - gyminus(2))
02708                                             / 2. / TOLERANCE;
02709                           dof_is_fixed[current_dof] = true;
02710                           current_dof++;
02711                           // Getting a 2nd order xyz is more tedious
02712                           Point nxmym = elem->point(n),
02713                                 nxmyp = elem->point(n),
02714                                 nxpym = elem->point(n),
02715                                 nxpyp = elem->point(n);
02716                           nxmym(0) -= TOLERANCE;
02717                           nxmym(1) -= TOLERANCE;
02718                           nxmyp(0) -= TOLERANCE;
02719                           nxmyp(1) += TOLERANCE;
02720                           nxpym(0) += TOLERANCE;
02721                           nxpym(1) -= TOLERANCE;
02722                           nxpyp(0) += TOLERANCE;
02723                           nxpyp(1) += TOLERANCE;
02724                           Gradient gxmym = g->component(var_component,
02725                                                         nxmym,
02726                                                         system.time);
02727                           Gradient gxmyp = g->component(var_component,
02728                                                         nxmyp,
02729                                                         system.time);
02730                           Gradient gxpym = g->component(var_component,
02731                                                         nxpym,
02732                                                         system.time);
02733                           Gradient gxpyp = g->component(var_component,
02734                                                         nxpyp,
02735                                                         system.time);
02736                           Number gxzplus = (gxpyp(2) - gxmyp(2))
02737                                          / 2. / TOLERANCE;
02738                           Number gxzminus = (gxpym(2) - gxmym(2))
02739                                           / 2. / TOLERANCE;
02740                           // xyz derivative
02741                           Ue(current_dof) = (gxzplus - gxzminus)
02742                                             / 2. / TOLERANCE;
02743                           dof_is_fixed[current_dof] = true;
02744                           current_dof++;
02745                         }
02746                     }
02747                 }
02748               // Assume that other C_ONE elements have a single nodal
02749               // value shape function and nodal gradient component
02750               // shape functions
02751               else if (cont == C_ONE)
02752                 {
02753                   libmesh_assert_equal_to (nc, 1 + dim);
02754                   Ue(current_dof) = f->component(var_component,
02755                                                  elem->point(n),
02756                                                  system.time);
02757                   dof_is_fixed[current_dof] = true;
02758                   current_dof++;
02759                   Gradient grad = g->component(var_component,
02760                                                elem->point(n),
02761                                                system.time);
02762                   for (unsigned int i=0; i!= dim; ++i)
02763                     {
02764                       Ue(current_dof) = grad(i);
02765                       dof_is_fixed[current_dof] = true;
02766                       current_dof++;
02767                     }
02768                 }
02769               else
02770                 libmesh_error();
02771             }
02772 
02773           // In 3D, project any edge values next
02774           if (dim > 2 && cont != DISCONTINUOUS)
02775             for (unsigned int e=0; e != elem->n_edges(); ++e)
02776               {
02777                 if (!is_boundary_edge[e])
02778                   continue;
02779 
02780                 FEInterface::dofs_on_edge(elem, dim, fe_type, e,
02781                                           side_dofs);
02782 
02783                 // Some edge dofs are on nodes and already
02784                 // fixed, others are free to calculate
02785                 unsigned int free_dofs = 0;
02786                 for (unsigned int i=0; i != side_dofs.size(); ++i)
02787                   if (!dof_is_fixed[side_dofs[i]])
02788                     free_dof[free_dofs++] = i;
02789 
02790                 // There may be nothing to project
02791                 if (!free_dofs)
02792                   continue;
02793 
02794                 Ke.resize (free_dofs, free_dofs); Ke.zero();
02795                 Fe.resize (free_dofs); Fe.zero();
02796                 // The new edge coefficients
02797                 DenseVector<Number> Uedge(free_dofs);
02798 
02799                 // Initialize FE data on the edge
02800                 fe->attach_quadrature_rule (qedgerule.get());
02801                 fe->edge_reinit (elem, e);
02802                 const unsigned int n_qp = qedgerule->n_points();
02803 
02804                 // Loop over the quadrature points
02805                 for (unsigned int qp=0; qp<n_qp; qp++)
02806                   {
02807                     // solution at the quadrature point
02808                     Number fineval = f->component(var_component,
02809                                                   xyz_values[qp],
02810                                                   system.time);
02811                     // solution grad at the quadrature point
02812                     Gradient finegrad;
02813                     if (cont == C_ONE)
02814                       finegrad = g->component(var_component,
02815                                               xyz_values[qp],
02816                                               system.time);
02817 
02818                     // Form edge projection matrix
02819                     for (unsigned int sidei=0, freei=0;
02820                          sidei != side_dofs.size(); ++sidei)
02821                       {
02822                         unsigned int i = side_dofs[sidei];
02823                         // fixed DoFs aren't test functions
02824                         if (dof_is_fixed[i])
02825                           continue;
02826                         for (unsigned int sidej=0, freej=0;
02827                              sidej != side_dofs.size(); ++sidej)
02828                           {
02829                             unsigned int j = side_dofs[sidej];
02830                             if (dof_is_fixed[j])
02831                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
02832                                            JxW[qp] * Ue(j);
02833                             else
02834                               Ke(freei,freej) += phi[i][qp] *
02835                                                  phi[j][qp] * JxW[qp];
02836                             if (cont == C_ONE)
02837                               {
02838                                 if (dof_is_fixed[j])
02839                                   Fe(freei) -= ((*dphi)[i][qp] *
02840                                                 (*dphi)[j][qp]) *
02841                                                 JxW[qp] * Ue(j);
02842                                 else
02843                                   Ke(freei,freej) += ((*dphi)[i][qp] *
02844                                                       (*dphi)[j][qp])
02845                                                       * JxW[qp];
02846                               }
02847                             if (!dof_is_fixed[j])
02848                               freej++;
02849                           }
02850                         Fe(freei) += phi[i][qp] * fineval * JxW[qp];
02851                         if (cont == C_ONE)
02852                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
02853                                        JxW[qp];
02854                         freei++;
02855                       }
02856                   }
02857 
02858                 Ke.cholesky_solve(Fe, Uedge);
02859 
02860                 // Transfer new edge solutions to element
02861                 for (unsigned int i=0; i != free_dofs; ++i)
02862                   {
02863                     Number &ui = Ue(side_dofs[free_dof[i]]);
02864                     libmesh_assert(std::abs(ui) < TOLERANCE ||
02865                            std::abs(ui - Uedge(i)) < TOLERANCE);
02866                     ui = Uedge(i);
02867                     dof_is_fixed[side_dofs[free_dof[i]]] = true;
02868                   }
02869               }
02870 
02871           // Project any side values (edges in 2D, faces in 3D)
02872           if (dim > 1 && cont != DISCONTINUOUS)
02873             for (unsigned int s=0; s != elem->n_sides(); ++s)
02874               {
02875                 if (!is_boundary_side[s])
02876                   continue;
02877 
02878                 FEInterface::dofs_on_side(elem, dim, fe_type, s,
02879                                           side_dofs);
02880 
02881                 // Some side dofs are on nodes/edges and already
02882                 // fixed, others are free to calculate
02883                 unsigned int free_dofs = 0;
02884                 for (unsigned int i=0; i != side_dofs.size(); ++i)
02885                   if (!dof_is_fixed[side_dofs[i]])
02886                     free_dof[free_dofs++] = i;
02887 
02888                 // There may be nothing to project
02889                 if (!free_dofs)
02890                   continue;
02891 
02892                 Ke.resize (free_dofs, free_dofs); Ke.zero();
02893                 Fe.resize (free_dofs); Fe.zero();
02894                 // The new side coefficients
02895                 DenseVector<Number> Uside(free_dofs);
02896 
02897                 // Initialize FE data on the side
02898                 fe->attach_quadrature_rule (qsiderule.get());
02899                 fe->reinit (elem, s);
02900                 const unsigned int n_qp = qsiderule->n_points();
02901 
02902                 // Loop over the quadrature points
02903                 for (unsigned int qp=0; qp<n_qp; qp++)
02904                   {
02905                     // solution at the quadrature point
02906                     Number fineval = f->component(var_component,
02907                                                   xyz_values[qp],
02908                                                   system.time);
02909                     // solution grad at the quadrature point
02910                     Gradient finegrad;
02911                     if (cont == C_ONE)
02912                       finegrad = g->component(var_component,
02913                                               xyz_values[qp],
02914                                               system.time);
02915 
02916                     // Form side projection matrix
02917                     for (unsigned int sidei=0, freei=0;
02918                          sidei != side_dofs.size(); ++sidei)
02919                       {
02920                         unsigned int i = side_dofs[sidei];
02921                         // fixed DoFs aren't test functions
02922                         if (dof_is_fixed[i])
02923                           continue;
02924                         for (unsigned int sidej=0, freej=0;
02925                              sidej != side_dofs.size(); ++sidej)
02926                           {
02927                             unsigned int j = side_dofs[sidej];
02928                             if (dof_is_fixed[j])
02929                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
02930                                            JxW[qp] * Ue(j);
02931                             else
02932                               Ke(freei,freej) += phi[i][qp] *
02933                                                  phi[j][qp] * JxW[qp];
02934                             if (cont == C_ONE)
02935                               {
02936                                 if (dof_is_fixed[j])
02937                                   Fe(freei) -= ((*dphi)[i][qp] *
02938                                                 (*dphi)[j][qp]) *
02939                                                JxW[qp] * Ue(j);
02940                                 else
02941                                   Ke(freei,freej) += ((*dphi)[i][qp] *
02942                                                       (*dphi)[j][qp])
02943                                                      * JxW[qp];
02944                               }
02945                             if (!dof_is_fixed[j])
02946                               freej++;
02947                           }
02948                         Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
02949                         if (cont == C_ONE)
02950                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
02951                                        JxW[qp];
02952                         freei++;
02953                       }
02954                   }
02955 
02956                 Ke.cholesky_solve(Fe, Uside);
02957 
02958                 // Transfer new side solutions to element
02959                 for (unsigned int i=0; i != free_dofs; ++i)
02960                   {
02961                     Number &ui = Ue(side_dofs[free_dof[i]]);
02962                     libmesh_assert(std::abs(ui) < TOLERANCE ||
02963                            std::abs(ui - Uside(i)) < TOLERANCE);
02964                     ui = Uside(i);
02965                     dof_is_fixed[side_dofs[free_dof[i]]] = true;
02966                   }
02967               }
02968 
02969           const dof_id_type
02970             first = new_vector.first_local_index(),
02971             last  = new_vector.last_local_index();
02972 
02973           // Lock the new_vector since it is shared among threads.
02974           {
02975             Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
02976 
02977             for (unsigned int i = 0; i < n_dofs; i++)
02978               if (dof_is_fixed[i] &&
02979                   (dof_indices[i] >= first) &&
02980                   (dof_indices[i] <  last))
02981                 {
02982                   new_vector.set(dof_indices[i], Ue(i));
02983                 }
02984           }
02985         }  // end elem loop
02986     } // end variables loop
02987 }


Member Data Documentation

Definition at line 203 of file system_projection.C.

Definition at line 205 of file system_projection.C.

Referenced by operator()().

const std::vector<unsigned int>& libMesh::BoundaryProjectSolution::variables [private]

Definition at line 204 of file system_projection.C.


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