// We now define the matrix assembly function for the // Biharmonic system. We need to first compute element // matrices and right-hand sides, and then take into // account the boundary conditions, which will be handled // via a penalty method. void assemble_biharmonic(EquationSystems& es, const std::string& system_name) { #ifdef LIBMESH_ENABLE_AMR #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES // It is a good idea to make sure we are assembling // the proper system. libmesh_assert_equal_to (system_name, "Biharmonic"); // Declare a performance log. Give it a descriptive // string to identify what part of the code we are // logging, since there may be many PerfLogs in an // application. PerfLog perf_log ("Matrix Assembly",false); // Get a constant reference to the mesh object. const MeshBase& mesh = es.get_mesh(); // The dimension that we are running const unsigned int dim = mesh.mesh_dimension(); // Get a reference to the LinearImplicitSystem we are solving LinearImplicitSystem& system = es.get_system<LinearImplicitSystem>("Biharmonic"); // A reference to the \p DofMap object for this system. The \p DofMap // object handles the index translation from node and element numbers // to degree of freedom numbers. We will talk more about the \p DofMap // in future examples. const DofMap& dof_map = system.get_dof_map(); // Get a constant reference to the Finite Element type // for the first (and only) variable in the system. FEType fe_type = dof_map.variable_type(0); // Build a Finite Element object of the specified type. Since the // \p FEBase::build() member dynamically creates memory we will // store the object as an \p UniquePtr<FEBase>. This can be thought // of as a pointer that will clean up after itself. UniquePtr<FEBase> fe (FEBase::build(dim, fe_type)); // Quadrature rule for numerical integration. // With 2D triangles, the Clough quadrature rule puts a Gaussian // quadrature rule on each of the 3 subelements UniquePtr<QBase> qrule(fe_type.default_quadrature_rule(dim)); // Tell the finite element object to use our quadrature rule. fe->attach_quadrature_rule (qrule.get()); // Declare a special finite element object for // boundary integration. UniquePtr<FEBase> fe_face (FEBase::build(dim, fe_type)); // Boundary integration requires another quadraure rule, // with dimensionality one less than the dimensionality // of the element. // In 1D, the Clough and Gauss quadrature rules are identical. UniquePtr<QBase> qface(fe_type.default_quadrature_rule(dim-1)); // Tell the finte element object to use our // quadrature rule. fe_face->attach_quadrature_rule (qface.get()); // Here we define some references to cell-specific data that // will be used to assemble the linear system. // We begin with the element Jacobian * quadrature weight at each // integration point. const std::vector<Real>& JxW = fe->get_JxW(); // The physical XY locations of the quadrature points on the element. // These might be useful for evaluating spatially varying material // properties at the quadrature points. const std::vector<Point>& q_point = fe->get_xyz(); // The element shape functions evaluated at the quadrature points. const std::vector<std::vector<Real> >& phi = fe->get_phi(); // The element shape function second derivatives evaluated at the // quadrature points. Note that for the simple biharmonic, shape // function first derivatives are unnecessary. const std::vector<std::vector<RealTensor> >& d2phi = fe->get_d2phi(); // For efficiency we will compute shape function laplacians n times, // not n^2 std::vector<Real> shape_laplacian; // Define data structures to contain the element matrix // and right-hand-side vector contribution. Following // basic finite element terminology we will denote these // "Ke" and "Fe". More detail is in example 3. DenseMatrix<Number> Ke; DenseVector<Number> Fe; // This vector will hold the degree of freedom indices for // the element. These define where in the global system // the element degrees of freedom get mapped. std::vector<dof_id_type> dof_indices; // Now we will loop over all the elements in the mesh. We will // compute the element matrix and right-hand-side contribution. See // example 3 for a discussion of the element iterators. MeshBase::const_element_iterator el = mesh.active_local_elements_begin(); const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end(); for ( ; el != end_el; ++el) { // Start logging the shape function initialization. // This is done through a simple function call with // the name of the event to log. perf_log.push("elem init"); // Store a pointer to the element we are currently // working on. This allows for nicer syntax later. const Elem* elem = *el; // Get the degree of freedom indices for the // current element. These define where in the global // matrix and right-hand-side this element will // contribute to. dof_map.dof_indices (elem, dof_indices); // Compute the element-specific data for the current // element. This involves computing the location of the // quadrature points (q_point) and the shape functions // (phi, dphi) for the current element. fe->reinit (elem); // Zero the element matrix and right-hand side before // summing them. Ke.resize (dof_indices.size(), dof_indices.size()); Fe.resize (dof_indices.size()); // Make sure there is enough room in this cache shape_laplacian.resize(dof_indices.size()); // Stop logging the shape function initialization. // If you forget to stop logging an event the PerfLog // object will probably catch the error and abort. perf_log.pop("elem init"); // Now we will build the element matrix. This involves // a double loop to integrate laplacians of the test funcions // (i) against laplacians of the trial functions (j). // // This step is why we need the Clough-Tocher elements - // these C1 differentiable elements have square-integrable // second derivatives. // // Now start logging the element matrix computation perf_log.push ("Ke"); for (unsigned int qp=0; qp<qrule->n_points(); qp++) { for (unsigned int i=0; i<phi.size(); i++) { shape_laplacian[i] = d2phi[i][qp](0,0)+d2phi[i][qp](1,1); if (dim == 3) shape_laplacian[i] += d2phi[i][qp](2,2); } for (unsigned int i=0; i<phi.size(); i++) for (unsigned int j=0; j<phi.size(); j++) Ke(i,j) += JxW[qp]* shape_laplacian[i]*shape_laplacian[j]; } // Stop logging the matrix computation perf_log.pop ("Ke"); // At this point the interior element integration has // been completed. However, we have not yet addressed // boundary conditions. For this example we will only // consider simple Dirichlet boundary conditions imposed // via the penalty method. Note that this is a fourth-order // problem: Dirichlet boundary conditions include *both* // boundary values and boundary normal fluxes. { // Start logging the boundary condition computation perf_log.push ("BCs"); // The penalty values, for solution boundary trace and flux. const Real penalty = 1e10; const Real penalty2 = 1e10; // The following loops over the sides of the element. // If the element has no neighbor on a side then that // side MUST live on a boundary of the domain. for (unsigned int s=0; s<elem->n_sides(); s++) if (elem->neighbor(s) == NULL) { // The value of the shape functions at the quadrature // points. const std::vector<std::vector<Real> >& phi_face = fe_face->get_phi(); // The value of the shape function derivatives at the // quadrature points. const std::vector<std::vector<RealGradient> >& dphi_face = fe_face->get_dphi(); // The Jacobian * Quadrature Weight at the quadrature // points on the face. const std::vector<Real>& JxW_face = fe_face->get_JxW(); // The XYZ locations (in physical space) of the // quadrature points on the face. This is where // we will interpolate the boundary value function. const std::vector<Point>& qface_point = fe_face->get_xyz(); const std::vector<Point>& face_normals = fe_face->get_normals(); // Compute the shape function values on the element // face. fe_face->reinit(elem, s); // Loop over the face quagrature points for integration. for (unsigned int qp=0; qp<qface->n_points(); qp++) { // The boundary value. Number value = exact_solution(qface_point[qp], es.parameters, "null", "void"); Gradient flux = exact_2D_derivative(qface_point[qp], es.parameters, "null", "void"); // Matrix contribution of the L2 projection. // Note that the basis function values are // integrated against test function values while // basis fluxes are integrated against test function // fluxes. for (unsigned int i=0; i<phi_face.size(); i++) for (unsigned int j=0; j<phi_face.size(); j++) Ke(i,j) += JxW_face[qp] * (penalty * phi_face[i][qp] * phi_face[j][qp] + penalty2 * (dphi_face[i][qp] * face_normals[qp]) * (dphi_face[j][qp] * face_normals[qp])); // Right-hand-side contribution of the L2 // projection. for (unsigned int i=0; i<phi_face.size(); i++) Fe(i) += JxW_face[qp] * (penalty * value * phi_face[i][qp] + penalty2 * (flux * face_normals[qp]) * (dphi_face[i][qp] * face_normals[qp])); } } // Stop logging the boundary condition computation perf_log.pop ("BCs"); } for (unsigned int qp=0; qp<qrule->n_points(); qp++) for (unsigned int i=0; i<phi.size(); i++) Fe(i) += JxW[qp]*phi[i][qp]*forcing_function(q_point[qp]); // The element matrix and right-hand-side are now built // for this element. Add them to the global matrix and // right-hand-side vector. The \p SparseMatrix::add_matrix() // and \p NumericVector::add_vector() members do this for us. // Start logging the insertion of the local (element) // matrix and vector into the global matrix and vector perf_log.push ("matrix insertion"); dof_map.constrain_element_matrix_and_vector(Ke, Fe, dof_indices); system.matrix->add_matrix (Ke, dof_indices); system.rhs->add_vector (Fe, dof_indices); // Stop logging the insertion of the local (element) // matrix and vector into the global matrix and vector perf_log.pop ("matrix insertion"); } // That's it. We don't need to do anything else to the // PerfLog. When it goes out of scope (at this function return) // it will print its log to the screen. Pretty easy, huh? #else #endif // #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES #endif // #ifdef LIBMESH_ENABLE_AMR }
void Biharmonic::JR::residual_and_jacobian(const NumericVector<Number> & u, NumericVector<Number> * R, SparseMatrix<Number> * J, NonlinearImplicitSystem &) { #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES if (!R && !J) return; // Declare a performance log. Give it a descriptive // string to identify what part of the code we are // logging, since there may be many PerfLogs in an // application. PerfLog perf_log ("Biharmonic Residual and Jacobian", false); // A reference to the DofMap object for this system. The DofMap // object handles the index translation from node and element numbers // to degree of freedom numbers. We will talk more about the DofMap // in future examples. const DofMap & dof_map = get_dof_map(); // Get a constant reference to the Finite Element type // for the first (and only) variable in the system. FEType fe_type = dof_map.variable_type(0); // Build a Finite Element object of the specified type. Since the // FEBase::build() member dynamically creates memory we will // store the object as a UniquePtr<FEBase>. This can be thought // of as a pointer that will clean up after itself. UniquePtr<FEBase> fe (FEBase::build(_biharmonic._dim, fe_type)); // Quadrature rule for numerical integration. // With 2D triangles, the Clough quadrature rule puts a Gaussian // quadrature rule on each of the 3 subelements UniquePtr<QBase> qrule(fe_type.default_quadrature_rule(_biharmonic._dim)); // Tell the finite element object to use our quadrature rule. fe->attach_quadrature_rule (qrule.get()); // Here we define some references to element-specific data that // will be used to assemble the linear system. // We begin with the element Jacobian * quadrature weight at each // integration point. const std::vector<Real> & JxW = fe->get_JxW(); // The element shape functions evaluated at the quadrature points. const std::vector<std::vector<Real> > & phi = fe->get_phi(); // The element shape functions' derivatives evaluated at the quadrature points. const std::vector<std::vector<RealGradient> > & dphi = fe->get_dphi(); // The element shape functions' second derivatives evaluated at the quadrature points. const std::vector<std::vector<RealTensor> > & d2phi = fe->get_d2phi(); // For efficiency we will compute shape function laplacians n times, // not n^2 std::vector<Real> Laplacian_phi_qp; // Define data structures to contain the element matrix // and right-hand-side vector contribution. Following // basic finite element terminology we will denote these // "Je" and "Re". More detail is in example 3. DenseMatrix<Number> Je; DenseVector<Number> Re; // This vector will hold the degree of freedom indices for // the element. These define where in the global system // the element degrees of freedom get mapped. std::vector<dof_id_type> dof_indices; // Old solution const NumericVector<Number> & u_old = *old_local_solution; // Now we will loop over all the elements in the mesh. We will // compute the element matrix and right-hand-side contribution. See // example 3 for a discussion of the element iterators. MeshBase::const_element_iterator el = _biharmonic._mesh.active_local_elements_begin(); const MeshBase::const_element_iterator end_el = _biharmonic._mesh.active_local_elements_end(); for ( ; el != end_el; ++el) { // Store a pointer to the element we are currently // working on. This allows for nicer syntax later. const Elem * elem = *el; // Get the degree of freedom indices for the // current element. These define where in the global // matrix and right-hand-side this element will // contribute to. dof_map.dof_indices (elem, dof_indices); // Compute the element-specific data for the current // element. This involves computing the location of the // quadrature points (q_point) and the shape function // values/derivatives (phi, dphi, d2phi) for the current element. fe->reinit (elem); // Zero the element matrix, the right-hand side and the Laplacian matrix // before summing them. if (J) Je.resize(dof_indices.size(), dof_indices.size()); if (R) Re.resize(dof_indices.size()); Laplacian_phi_qp.resize(dof_indices.size()); for (unsigned int qp=0; qp<qrule->n_points(); qp++) { // AUXILIARY QUANTITIES: // Residual and Jacobian share a few calculations: // at the very least, in the case of interfacial energy only with a constant mobility, // both calculations use Laplacian_phi_qp; more is shared the case of a concentration-dependent // mobility and bulk potentials. Number u_qp = 0.0, u_old_qp = 0.0, Laplacian_u_qp = 0.0, Laplacian_u_old_qp = 0.0; Gradient grad_u_qp(0.0, 0.0, 0.0), grad_u_old_qp(0.0, 0.0, 0.0); Number M_qp = 1.0, M_old_qp = 1.0, M_prime_qp = 0.0, M_prime_old_qp = 0.0; for (unsigned int i=0; i<phi.size(); i++) { Laplacian_phi_qp[i] = d2phi[i][qp](0, 0); grad_u_qp(0) += u(dof_indices[i])*dphi[i][qp](0); grad_u_old_qp(0) += u_old(dof_indices[i])*dphi[i][qp](0); if (_biharmonic._dim > 1) { Laplacian_phi_qp[i] += d2phi[i][qp](1, 1); grad_u_qp(1) += u(dof_indices[i])*dphi[i][qp](1); grad_u_old_qp(1) += u_old(dof_indices[i])*dphi[i][qp](1); } if (_biharmonic._dim > 2) { Laplacian_phi_qp[i] += d2phi[i][qp](2, 2); grad_u_qp(2) += u(dof_indices[i])*dphi[i][qp](2); grad_u_old_qp(2) += u_old(dof_indices[i])*dphi[i][qp](2); } u_qp += phi[i][qp]*u(dof_indices[i]); u_old_qp += phi[i][qp]*u_old(dof_indices[i]); Laplacian_u_qp += Laplacian_phi_qp[i]*u(dof_indices[i]); Laplacian_u_old_qp += Laplacian_phi_qp[i]*u_old(dof_indices[i]); } // for i if (_biharmonic._degenerate) { M_qp = 1.0 - u_qp*u_qp; M_old_qp = 1.0 - u_old_qp*u_old_qp; M_prime_qp = -2.0*u_qp; M_prime_old_qp = -2.0*u_old_qp; } // ELEMENT RESIDUAL AND JACOBIAN for (unsigned int i=0; i<phi.size(); i++) { // RESIDUAL if (R) { Number ri = 0.0, ri_old = 0.0; ri -= Laplacian_phi_qp[i]*M_qp*_biharmonic._kappa*Laplacian_u_qp; ri_old -= Laplacian_phi_qp[i]*M_old_qp*_biharmonic._kappa*Laplacian_u_old_qp; if (_biharmonic._degenerate) { ri -= (dphi[i][qp]*grad_u_qp)*M_prime_qp*(_biharmonic._kappa*Laplacian_u_qp); ri_old -= (dphi[i][qp]*grad_u_old_qp)*M_prime_old_qp*(_biharmonic._kappa*Laplacian_u_old_qp); } if (_biharmonic._cahn_hillard) { if (_biharmonic._energy == DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_WELL) { ri += Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp; ri_old += Laplacian_phi_qp[i]*M_old_qp*_biharmonic._theta_c*(u_old_qp*u_old_qp - 1.0)*u_old_qp; if (_biharmonic._degenerate) { ri += (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp; ri_old += (dphi[i][qp]*grad_u_old_qp)*M_prime_old_qp*_biharmonic._theta_c*(u_old_qp*u_old_qp - 1.0)*u_old_qp; } }// if (_biharmonic._energy == DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_WELL) if (_biharmonic._energy == DOUBLE_OBSTACLE || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) { ri -= Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*u_qp; ri_old -= Laplacian_phi_qp[i]*M_old_qp*_biharmonic._theta_c*u_old_qp; if (_biharmonic._degenerate) { ri -= (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*u_qp; ri_old -= (dphi[i][qp]*grad_u_old_qp)*M_prime_old_qp*_biharmonic._theta_c*u_old_qp; } } // if (_biharmonic._energy == DOUBLE_OBSTACLE || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) if (_biharmonic._energy == LOG_DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) { switch(_biharmonic._log_truncation) { case 2: break; case 3: break; default: break; }// switch(_biharmonic._log_truncation) }// if (_biharmonic._energy == LOG_DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) }// if (_biharmonic._cahn_hillard) Re(i) += JxW[qp]*((u_qp-u_old_qp)*phi[i][qp]-_biharmonic._dt*0.5*((2.0-_biharmonic._cnWeight)*ri + _biharmonic._cnWeight*ri_old)); } // if (R) // JACOBIAN if (J) { Number M_prime_prime_qp = 0.0; if (_biharmonic._degenerate) M_prime_prime_qp = -2.0; for (unsigned int j=0; j<phi.size(); j++) { Number ri_j = 0.0; ri_j -= Laplacian_phi_qp[i]*M_qp*_biharmonic._kappa*Laplacian_phi_qp[j]; if (_biharmonic._degenerate) { ri_j -= Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._kappa*Laplacian_u_qp + (dphi[i][qp]*dphi[j][qp])*M_prime_qp*(_biharmonic._kappa*Laplacian_u_qp) + (dphi[i][qp]*grad_u_qp)*(M_prime_prime_qp*phi[j][qp])*(_biharmonic._kappa*Laplacian_u_qp) + (dphi[i][qp]*grad_u_qp)*(M_prime_qp)*(_biharmonic._kappa*Laplacian_phi_qp[j]); } if (_biharmonic._cahn_hillard) { if (_biharmonic._energy == DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_WELL) { ri_j += Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp + Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*(3.0*u_qp*u_qp - 1.0)*phi[j][qp] + (dphi[i][qp]*dphi[j][qp])*M_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp + (dphi[i][qp]*grad_u_qp)*M_prime_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp + (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*(3.0*u_qp*u_qp - 1.0)*phi[j][qp]; }// if (_biharmonic._energy == DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_WELL) if (_biharmonic._energy == DOUBLE_OBSTACLE || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) { ri_j -= Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._theta_c*u_qp + Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*phi[j][qp] + (dphi[i][qp]*dphi[j][qp])*M_prime_qp*_biharmonic._theta_c*u_qp + (dphi[i][qp]*grad_u_qp)*M_prime_prime_qp*_biharmonic._theta_c*u_qp + (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*phi[j][qp]; } // if (_biharmonic._energy == DOUBLE_OBSTACLE || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) if (_biharmonic._energy == LOG_DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) { switch(_biharmonic._log_truncation) { case 2: break; case 3: break; default: break; }// switch(_biharmonic._log_truncation) }// if (_biharmonic._energy == LOG_DOUBLE_WELL || _biharmonic._energy == LOG_DOUBLE_OBSTACLE) }// if (_biharmonic._cahn_hillard) Je(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] - 0.5*_biharmonic._dt*(2.0-_biharmonic._cnWeight)*ri_j); } // for j } // if (J) } // for i } // for qp // The element matrix and right-hand-side are now built // for this element. Add them to the global matrix and // right-hand-side vector. The SparseMatrix::add_matrix() // and NumericVector::add_vector() members do this for us. // Start logging the insertion of the local (element) // matrix and vector into the global matrix and vector if (R) { // If the mesh has hanging nodes (e.g., as a result of refinement), those need to be constrained. dof_map.constrain_element_vector(Re, dof_indices); R->add_vector(Re, dof_indices); } if (J) { // If the mesh has hanging nodes (e.g., as a result of refinement), those need to be constrained. dof_map.constrain_element_matrix(Je, dof_indices); J->add_matrix(Je, dof_indices); } } // for el #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES }
void ExactErrorEstimator::estimate_error (const System & system, ErrorVector & error_per_cell, const NumericVector<Number> * solution_vector, bool estimate_parent_error) { // Ignore the fact that this variable is unused when !LIBMESH_ENABLE_AMR libmesh_ignore(estimate_parent_error); // The current mesh const MeshBase & mesh = system.get_mesh(); // The dimensionality of the mesh const unsigned int dim = mesh.mesh_dimension(); // The number of variables in the system const unsigned int n_vars = system.n_vars(); // The DofMap for this system const DofMap & dof_map = system.get_dof_map(); // Resize the error_per_cell vector to be // the number of elements, initialize it to 0. error_per_cell.resize (mesh.max_elem_id()); std::fill (error_per_cell.begin(), error_per_cell.end(), 0.); // Prepare current_local_solution to localize a non-standard // solution vector if necessary if (solution_vector && solution_vector != system.solution.get()) { NumericVector<Number> * newsol = const_cast<NumericVector<Number> *>(solution_vector); System & sys = const_cast<System &>(system); newsol->swap(*sys.solution); sys.update(); } // Loop over all the variables in the system for (unsigned int var=0; var<n_vars; var++) { // Possibly skip this variable if (error_norm.weight(var) == 0.0) continue; // The (string) name of this variable const std::string & var_name = system.variable_name(var); // The type of finite element to use for this variable const FEType & fe_type = dof_map.variable_type (var); UniquePtr<FEBase> fe (FEBase::build (dim, fe_type)); // Build an appropriate Gaussian quadrature rule UniquePtr<QBase> qrule = fe_type.default_quadrature_rule (dim, _extra_order); fe->attach_quadrature_rule (qrule.get()); // Prepare a global solution and a MeshFunction of the fine system if we need one UniquePtr<MeshFunction> fine_values; UniquePtr<NumericVector<Number> > fine_soln = NumericVector<Number>::build(system.comm()); if (_equation_systems_fine) { const System & fine_system = _equation_systems_fine->get_system(system.name()); std::vector<Number> global_soln; // FIXME - we're assuming that the fine system solution gets // used even when a different vector is used for the coarse // system fine_system.update_global_solution(global_soln); fine_soln->init (cast_int<numeric_index_type>(global_soln.size()), true, SERIAL); (*fine_soln) = global_soln; fine_values = UniquePtr<MeshFunction> (new MeshFunction(*_equation_systems_fine, *fine_soln, fine_system.get_dof_map(), fine_system.variable_number(var_name))); fine_values->init(); } else { // Initialize functors if we're using them for (unsigned int i=0; i != _exact_values.size(); ++i) if (_exact_values[i]) _exact_values[i]->init(); for (unsigned int i=0; i != _exact_derivs.size(); ++i) if (_exact_derivs[i]) _exact_derivs[i]->init(); for (unsigned int i=0; i != _exact_hessians.size(); ++i) if (_exact_hessians[i]) _exact_hessians[i]->init(); } // Request the data we'll need to compute with fe->get_JxW(); fe->get_phi(); fe->get_dphi(); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES fe->get_d2phi(); #endif fe->get_xyz(); #ifdef LIBMESH_ENABLE_AMR // If we compute on parent elements, we'll want to do so only // once on each, so we need to keep track of which we've done. std::vector<bool> computed_var_on_parent; if (estimate_parent_error) computed_var_on_parent.resize(error_per_cell.size(), false); #endif // TODO: this ought to be threaded (and using subordinate // MeshFunction objects in each thread rather than a single // master) // Iterate over all the active elements in the mesh // that live on this processor. MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); for (;elem_it != elem_end; ++elem_it) { // e is necessarily an active element on the local processor const Elem * elem = *elem_it; const dof_id_type e_id = elem->id(); #ifdef LIBMESH_ENABLE_AMR // See if the parent of element e has been examined yet; // if not, we may want to compute the estimator on it const Elem * parent = elem->parent(); // We only can compute and only need to compute on // parents with all active children bool compute_on_parent = true; if (!parent || !estimate_parent_error) compute_on_parent = false; else for (unsigned int c=0; c != parent->n_children(); ++c) if (!parent->child_ptr(c)->active()) compute_on_parent = false; if (compute_on_parent && !computed_var_on_parent[parent->id()]) { computed_var_on_parent[parent->id()] = true; // Compute a projection onto the parent DenseVector<Number> Uparent; FEBase::coarsened_dof_values(*(system.current_local_solution), dof_map, parent, Uparent, var, false); error_per_cell[parent->id()] += static_cast<ErrorVectorReal> (find_squared_element_error(system, var_name, parent, Uparent, fe.get(), fine_values.get())); } #endif // Get the local to global degree of freedom maps std::vector<dof_id_type> dof_indices; dof_map.dof_indices (elem, dof_indices, var); const unsigned int n_dofs = cast_int<unsigned int>(dof_indices.size()); DenseVector<Number> Uelem(n_dofs); for (unsigned int i=0; i != n_dofs; ++i) Uelem(i) = system.current_solution(dof_indices[i]); error_per_cell[e_id] += static_cast<ErrorVectorReal> (find_squared_element_error(system, var_name, elem, Uelem, fe.get(), fine_values.get())); } // End loop over active local elements } // End loop over variables // Each processor has now computed the error contribuions // for its local elements. We need to sum the vector // and then take the square-root of each component. Note // that we only need to sum if we are running on multiple // processors, and we only need to take the square-root // if the value is nonzero. There will in general be many // zeros for the inactive elements. // First sum the vector of estimated error values this->reduce_error(error_per_cell, system.comm()); // Compute the square-root of each component. { LOG_SCOPE("std::sqrt()", "ExactErrorEstimator"); for (dof_id_type i=0; i<error_per_cell.size(); i++) if (error_per_cell[i] != 0.) { libmesh_assert_greater (error_per_cell[i], 0.); error_per_cell[i] = std::sqrt(error_per_cell[i]); } } // If we used a non-standard solution before, now is the time to fix // the current_local_solution if (solution_vector && solution_vector != system.solution.get()) { NumericVector<Number> * newsol = const_cast<NumericVector<Number> *>(solution_vector); System & sys = const_cast<System &>(system); newsol->swap(*sys.solution); sys.update(); } }
// We now define the matrix and rhs vector assembly function // for the shell system. This function implements the // linear Kirchhoff-Love theory for thin shells. At the // end we also take into account the boundary conditions // here, using the penalty method. void assemble_shell (EquationSystems& es, const std::string& system_name) { // It is a good idea to make sure we are assembling // the proper system. libmesh_assert_equal_to (system_name, "Shell"); // Get a constant reference to the mesh object. const MeshBase& mesh = es.get_mesh(); // Get a reference to the shell system object. LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Shell"); // Get the shell parameters that we need during assembly. const Real h = es.parameters.get<Real> ("thickness"); const Real E = es.parameters.get<Real> ("young's modulus"); const Real nu = es.parameters.get<Real> ("poisson ratio"); const Real q = es.parameters.get<Real> ("uniform load"); // Compute the membrane stiffness \p K and the bending // rigidity \p D from these parameters. const Real K = E * h / (1-nu*nu); const Real D = E * h*h*h / (12*(1-nu*nu)); // Numeric ids corresponding to each variable in the system. const unsigned int u_var = system.variable_number ("u"); const unsigned int v_var = system.variable_number ("v"); const unsigned int w_var = system.variable_number ("w"); // Get the Finite Element type for "u". Note this will be // the same as the type for "v" and "w". FEType fe_type = system.variable_type (u_var); // Build a Finite Element object of the specified type. UniquePtr<FEBase> fe (FEBase::build(2, fe_type)); // A Gauss quadrature rule for numerical integration. // For subdivision shell elements, a single Gauss point per // element is sufficient, hence we use extraorder = 0. const int extraorder = 0; UniquePtr<QBase> qrule (fe_type.default_quadrature_rule (2, extraorder)); // Tell the finite element object to use our quadrature rule. fe->attach_quadrature_rule (qrule.get()); // The element Jacobian * quadrature weight at each integration point. const std::vector<Real>& JxW = fe->get_JxW(); // The surface tangents in both directions at the quadrature points. const std::vector<RealGradient>& dxyzdxi = fe->get_dxyzdxi(); const std::vector<RealGradient>& dxyzdeta = fe->get_dxyzdeta(); // The second partial derivatives at the quadrature points. const std::vector<RealGradient>& d2xyzdxi2 = fe->get_d2xyzdxi2(); const std::vector<RealGradient>& d2xyzdeta2 = fe->get_d2xyzdeta2(); const std::vector<RealGradient>& d2xyzdxideta = fe->get_d2xyzdxideta(); // The element shape function and its derivatives evaluated at the // quadrature points. const std::vector<std::vector<Real> >& phi = fe->get_phi(); const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi(); const std::vector<std::vector<RealTensor> >& d2phi = fe->get_d2phi(); // A reference to the \p DofMap object for this system. The \p DofMap // object handles the index translation from node and element numbers // to degree of freedom numbers. const DofMap & dof_map = system.get_dof_map(); // Define data structures to contain the element stiffness matrix // and right-hand-side vector contribution. Following // basic finite element terminology we will denote these // "Ke" and "Fe". DenseMatrix<Number> Ke; DenseVector<Number> Fe; DenseSubMatrix<Number> Kuu(Ke), Kuv(Ke), Kuw(Ke), Kvu(Ke), Kvv(Ke), Kvw(Ke), Kwu(Ke), Kwv(Ke), Kww(Ke); DenseSubVector<Number> Fu(Fe), Fv(Fe), Fw(Fe); // This vector will hold the degree of freedom indices for // the element. These define where in the global system // the element degrees of freedom get mapped. std::vector<dof_id_type> dof_indices; std::vector<dof_id_type> dof_indices_u; std::vector<dof_id_type> dof_indices_v; std::vector<dof_id_type> dof_indices_w; // Now we will loop over all the elements in the mesh. We will // compute the element matrix and right-hand-side contribution. MeshBase::const_element_iterator el = mesh.active_local_elements_begin(); const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end(); for (; el != end_el; ++el) { // Store a pointer to the element we are currently // working on. This allows for nicer syntax later. const Elem* elem = *el; // The ghost elements at the boundaries need to be excluded // here, as they don't belong to the physical shell, // but serve for a proper boundary treatment only. libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION); const Tri3Subdivision* sd_elem = static_cast<const Tri3Subdivision*> (elem); if (sd_elem->is_ghost()) continue; // Get the degree of freedom indices for the // current element. These define where in the global // matrix and right-hand-side this element will // contribute to. dof_map.dof_indices (elem, dof_indices); dof_map.dof_indices (elem, dof_indices_u, u_var); dof_map.dof_indices (elem, dof_indices_v, v_var); dof_map.dof_indices (elem, dof_indices_w, w_var); const std::size_t n_dofs = dof_indices.size(); const std::size_t n_u_dofs = dof_indices_u.size(); const std::size_t n_v_dofs = dof_indices_v.size(); const std::size_t n_w_dofs = dof_indices_w.size(); // Compute the element-specific data for the current // element. This involves computing the location of the // quadrature points and the shape functions // (phi, dphi, d2phi) for the current element. fe->reinit (elem); // Zero the element matrix and right-hand side before // summing them. We use the resize member here because // the number of degrees of freedom might have changed from // the last element. Ke.resize (n_dofs, n_dofs); Fe.resize (n_dofs); // Reposition the submatrices... The idea is this: // // - - - - // | Kuu Kuv Kuw | | Fu | // Ke = | Kvu Kvv Kvw |; Fe = | Fv | // | Kwu Kwv Kww | | Fw | // - - - - // // The \p DenseSubMatrix.repostition () member takes the // (row_offset, column_offset, row_size, column_size). // // Similarly, the \p DenseSubVector.reposition () member // takes the (row_offset, row_size) Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs); Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs); Kuw.reposition (u_var*n_u_dofs, w_var*n_u_dofs, n_u_dofs, n_w_dofs); Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs); Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs); Kvw.reposition (v_var*n_v_dofs, w_var*n_v_dofs, n_v_dofs, n_w_dofs); Kwu.reposition (w_var*n_w_dofs, u_var*n_w_dofs, n_w_dofs, n_u_dofs); Kwv.reposition (w_var*n_w_dofs, v_var*n_w_dofs, n_w_dofs, n_v_dofs); Kww.reposition (w_var*n_w_dofs, w_var*n_w_dofs, n_w_dofs, n_w_dofs); Fu.reposition (u_var*n_u_dofs, n_u_dofs); Fv.reposition (v_var*n_u_dofs, n_v_dofs); Fw.reposition (w_var*n_u_dofs, n_w_dofs); // Now we will build the element matrix and right-hand-side. for (unsigned int qp=0; qp<qrule->n_points(); ++qp) { // First, we compute the external force resulting // from a load q distributed uniformly across the plate. // Since the load is supposed to be transverse to the plate, // it affects the z-direction, i.e. the "w" variable. for (unsigned int i=0; i<n_u_dofs; ++i) Fw(i) += JxW[qp] * phi[i][qp] * q; // Next, we assemble the stiffness matrix. This is only valid // for the linear theory, i.e., for small deformations, where // reference and deformed surface metrics are indistinguishable. // Get the three surface basis vectors. const RealVectorValue & a1 = dxyzdxi[qp]; const RealVectorValue & a2 = dxyzdeta[qp]; RealVectorValue a3 = a1.cross(a2); const Real jac = a3.size(); // the surface Jacobian libmesh_assert_greater (jac, 0); a3 /= jac; // the shell director a3 is normalized to unit length // Get the derivatives of the surface tangents. const RealVectorValue & a11 = d2xyzdxi2[qp]; const RealVectorValue & a22 = d2xyzdeta2[qp]; const RealVectorValue & a12 = d2xyzdxideta[qp]; // Compute the three covariant components of the first // fundamental form of the surface. const RealVectorValue a(a1*a1, a2*a2, a1*a2); // The elastic H matrix in Voigt's notation, computed from the // covariant components of the first fundamental form rather // than the contravariant components, exploiting that the // contravariant first fundamental form is the inverse of the // covatiant first fundamental form (hence the determinant etc.). RealTensorValue H; H(0,0) = a(1) * a(1); H(0,1) = H(1,0) = nu * a(1) * a(0) + (1-nu) * a(2) * a(2); H(0,2) = H(2,0) = -a(1) * a(2); H(1,1) = a(0) * a(0); H(1,2) = H(2,1) = -a(0) * a(2); H(2,2) = 0.5 * ((1-nu) * a(1) * a(0) + (1+nu) * a(2) * a(2)); const Real det = a(0) * a(1) - a(2) * a(2); libmesh_assert_not_equal_to (det * det, 0); H /= det * det; // Precompute come cross products for the bending part below. const RealVectorValue a11xa2 = a11.cross(a2); const RealVectorValue a22xa2 = a22.cross(a2); const RealVectorValue a12xa2 = a12.cross(a2); const RealVectorValue a1xa11 = a1.cross(a11); const RealVectorValue a1xa22 = a1.cross(a22); const RealVectorValue a1xa12 = a1.cross(a12); const RealVectorValue a2xa3 = a2.cross(a3); const RealVectorValue a3xa1 = a3.cross(a1); // Loop over all pairs of nodes I,J. for (unsigned int i=0; i<n_u_dofs; ++i) { for (unsigned int j=0; j<n_u_dofs; ++j) { // The membrane strain matrices in Voigt's notation. RealTensorValue MI, MJ; for (unsigned int k=0; k<3; ++k) { MI(0,k) = dphi[i][qp](0) * a1(k); MI(1,k) = dphi[i][qp](1) * a2(k); MI(2,k) = dphi[i][qp](1) * a1(k) + dphi[i][qp](0) * a2(k); MJ(0,k) = dphi[j][qp](0) * a1(k); MJ(1,k) = dphi[j][qp](1) * a2(k); MJ(2,k) = dphi[j][qp](1) * a1(k) + dphi[j][qp](0) * a2(k); } // The bending strain matrices in Voigt's notation. RealTensorValue BI, BJ; for (unsigned int k=0; k<3; ++k) { const Real term_ik = dphi[i][qp](0) * a2xa3(k) + dphi[i][qp](1) * a3xa1(k); BI(0,k) = -d2phi[i][qp](0,0) * a3(k) +(dphi[i][qp](0) * a11xa2(k) + dphi[i][qp](1) * a1xa11(k) + (a3*a11) * term_ik) / jac; BI(1,k) = -d2phi[i][qp](1,1) * a3(k) +(dphi[i][qp](0) * a22xa2(k) + dphi[i][qp](1) * a1xa22(k) + (a3*a22) * term_ik) / jac; BI(2,k) = 2 * (-d2phi[i][qp](0,1) * a3(k) +(dphi[i][qp](0) * a12xa2(k) + dphi[i][qp](1) * a1xa12(k) + (a3*a12) * term_ik) / jac); const Real term_jk = dphi[j][qp](0) * a2xa3(k) + dphi[j][qp](1) * a3xa1(k); BJ(0,k) = -d2phi[j][qp](0,0) * a3(k) +(dphi[j][qp](0) * a11xa2(k) + dphi[j][qp](1) * a1xa11(k) + (a3*a11) * term_jk) / jac; BJ(1,k) = -d2phi[j][qp](1,1) * a3(k) +(dphi[j][qp](0) * a22xa2(k) + dphi[j][qp](1) * a1xa22(k) + (a3*a22) * term_jk) / jac; BJ(2,k) = 2 * (-d2phi[j][qp](0,1) * a3(k) +(dphi[j][qp](0) * a12xa2(k) + dphi[j][qp](1) * a1xa12(k) + (a3*a12) * term_jk) / jac); } // The total stiffness matrix coupling the nodes // I and J is a sum of membrane and bending // contributions according to the following formula. const RealTensorValue KIJ = JxW[qp] * K * MI.transpose() * H * MJ + JxW[qp] * D * BI.transpose() * H * BJ; // Insert the components of the coupling stiffness // matrix \p KIJ into the corresponding directional // submatrices. Kuu(i,j) += KIJ(0,0); Kuv(i,j) += KIJ(0,1); Kuw(i,j) += KIJ(0,2); Kvu(i,j) += KIJ(1,0); Kvv(i,j) += KIJ(1,1); Kvw(i,j) += KIJ(1,2); Kwu(i,j) += KIJ(2,0); Kwv(i,j) += KIJ(2,1); Kww(i,j) += KIJ(2,2); } } } // end of the quadrature point qp-loop // The element matrix and right-hand-side are now built // for this element. Add them to the global matrix and // right-hand-side vector. The \p NumericMatrix::add_matrix() // and \p NumericVector::add_vector() members do this for us. system.matrix->add_matrix (Ke, dof_indices); system.rhs->add_vector (Fe, dof_indices); } // end of non-ghost element loop // Next, we apply the boundary conditions. In this case, // all boundaries are clamped by the penalty method, using // the special "ghost" nodes along the boundaries. Note // that there are better ways to implement boundary conditions // for subdivision shells. We use the simplest way here, // which is known to be overly restrictive and will lead to // a slightly too small deformation of the plate. el = mesh.active_local_elements_begin(); for (; el != end_el; ++el) { // Store a pointer to the element we are currently // working on. This allows for nicer syntax later. const Elem* elem = *el; // For the boundary conditions, we only need to loop over // the ghost elements. libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION); const Tri3Subdivision* gh_elem = static_cast<const Tri3Subdivision*> (elem); if (!gh_elem->is_ghost()) continue; // Find the side which is part of the physical plate boundary, // that is, the boundary of the original mesh without ghosts. for (unsigned int s=0; s<elem->n_sides(); ++s) { const Tri3Subdivision* nb_elem = static_cast<const Tri3Subdivision*> (elem->neighbor(s)); if (nb_elem == NULL || nb_elem->is_ghost()) continue; /* * Determine the four nodes involved in the boundary * condition treatment of this side. The \p MeshTools::Subdiv * namespace provides lookup tables \p next and \p prev * for an efficient determination of the next and previous * nodes of an element, respectively. * * n4 * / \ * / gh \ * n2 ---- n3 * \ nb / * \ / * n1 */ Node* nodes [4]; // n1, n2, n3, n4 nodes[1] = gh_elem->get_node(s); // n2 nodes[2] = gh_elem->get_node(MeshTools::Subdivision::next[s]); // n3 nodes[3] = gh_elem->get_node(MeshTools::Subdivision::prev[s]); // n4 // The node in the interior of the domain, \p n1, is the // hardest to find. Walk along the edges of element \p nb until // we have identified it. unsigned int n_int = 0; nodes[0] = nb_elem->get_node(0); while (nodes[0]->id() == nodes[1]->id() || nodes[0]->id() == nodes[2]->id()) nodes[0] = nb_elem->get_node(++n_int); // The penalty value. \f$ \frac{1}{\epsilon} \f$ const Real penalty = 1.e10; // With this simple method, clamped boundary conditions are // obtained by penalizing the displacements of all four nodes. // This ensures that the displacement field vanishes on the // boundary side \p s. for (unsigned int n=0; n<4; ++n) { const dof_id_type u_dof = nodes[n]->dof_number (system.number(), u_var, 0); const dof_id_type v_dof = nodes[n]->dof_number (system.number(), v_var, 0); const dof_id_type w_dof = nodes[n]->dof_number (system.number(), w_var, 0); system.matrix->add (u_dof, u_dof, penalty); system.matrix->add (v_dof, v_dof, penalty); system.matrix->add (w_dof, w_dof, penalty); } } } // end of ghost element loop }