2 #include "libmesh/mesh.h" 
    3 #include "libmesh/quadrature.h" 
    4 #include "libmesh/dense_matrix.h" 
    5 #include "libmesh/dense_vector.h" 
    6 #include "libmesh/sparse_matrix.h" 
    7 #include "libmesh/fourth_error_estimators.h" 
    8 #include "libmesh/dof_map.h" 
    9 #include "libmesh/numeric_vector.h" 
   10 #include "libmesh/periodic_boundaries.h" 
   11 #include "libmesh/periodic_boundary.h" 
   12 #include "libmesh/elem.h" 
   20                    const std::string & name_in,
 
   21                    const unsigned int number_in) :
 
   26 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES 
   27   libmesh_error_msg(
"Must have second derivatives enabled");
 
   30 #ifdef LIBMESH_ENABLE_PERIODIC 
   32   DofMap & dof_map = get_dof_map();
 
   74 #endif // LIBMESH_ENABLE_PERIODIC 
   81   attach_init_object(*
this);
 
   84   nonlinear_solver->residual_and_jacobian_object = 
this;
 
   87   nonlinear_solver->bounds_object = 
this;
 
   96   if (_biharmonic._verbose)
 
  113   *(old_local_solution) = *(current_local_solution);
 
  115   if (_biharmonic._verbose)
 
  134   return (r < width) ? 1.0 : -0.5;
 
  153   return (r < width) ? 1.0 : -0.5;
 
  168   Real r = 
sqrt((p(0)-center(0))*(p(0)-center(0)));
 
  169   return (r < width) ? 1.0 : -0.5;
 
  193 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  201   PerfLog perf_log (
"Biharmonic Residual and Jacobian", 
false);
 
  207   const DofMap & dof_map = get_dof_map();
 
  217   std::unique_ptr<FEBase> fe (FEBase::build(_biharmonic._dim, fe_type));
 
  225   fe->attach_quadrature_rule (qrule.get());
 
  231   const std::vector<Real> & JxW = fe->get_JxW();
 
  234   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  237   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  240   const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
 
  244   std::vector<Real> Laplacian_phi_qp;
 
  256   std::vector<dof_id_type> dof_indices;
 
  264   for (
const auto & elem : _biharmonic._mesh.active_local_element_ptr_range())
 
  272       const unsigned int n_dofs =
 
  273         cast_int<unsigned int>(dof_indices.size());
 
  289       Laplacian_phi_qp.resize(
n_dofs);
 
  291       for (
unsigned int qp=0; qp<qrule->n_points(); qp++)
 
  301             Laplacian_u_qp = 0.0,
 
  302             Laplacian_u_old_qp = 0.0;
 
  305             grad_u_qp(0.0, 0.0, 0.0),
 
  306             grad_u_old_qp(0.0, 0.0, 0.0);
 
  312             M_prime_old_qp = 0.0;
 
  314           for (
unsigned int i=0; i<
n_dofs; i++)
 
  316               Laplacian_phi_qp[i] = d2phi[i][qp](0, 0);
 
  317               grad_u_qp(0) += u(dof_indices[i])*dphi[i][qp](0);
 
  318               grad_u_old_qp(0) += u_old(dof_indices[i])*dphi[i][qp](0);
 
  320               if (_biharmonic._dim > 1)
 
  322                   Laplacian_phi_qp[i] += d2phi[i][qp](1, 1);
 
  323                   grad_u_qp(1) += u(dof_indices[i])*dphi[i][qp](1);
 
  324                   grad_u_old_qp(1) += u_old(dof_indices[i])*dphi[i][qp](1);
 
  326               if (_biharmonic._dim > 2)
 
  328                   Laplacian_phi_qp[i] += d2phi[i][qp](2, 2);
 
  329                   grad_u_qp(2) += u(dof_indices[i])*dphi[i][qp](2);
 
  330                   grad_u_old_qp(2) += u_old(dof_indices[i])*dphi[i][qp](2);
 
  332               u_qp     += phi[i][qp]*u(dof_indices[i]);
 
  333               u_old_qp += phi[i][qp]*u_old(dof_indices[i]);
 
  334               Laplacian_u_qp     += Laplacian_phi_qp[i]*u(dof_indices[i]);
 
  335               Laplacian_u_old_qp += Laplacian_phi_qp[i]*u_old(dof_indices[i]);
 
  338           if (_biharmonic._degenerate)
 
  340               M_qp           = 1.0 - u_qp*u_qp;
 
  341               M_old_qp       = 1.0 - u_old_qp*u_old_qp;
 
  342               M_prime_qp     = -2.0*u_qp;
 
  343               M_prime_old_qp = -2.0*u_old_qp;
 
  347           for (
unsigned int i=0; i<
n_dofs; i++)
 
  352                   Number ri = 0.0, ri_old = 0.0;
 
  353                   ri     -= Laplacian_phi_qp[i]*M_qp*_biharmonic._kappa*Laplacian_u_qp;
 
  354                   ri_old -= Laplacian_phi_qp[i]*M_old_qp*_biharmonic._kappa*Laplacian_u_old_qp;
 
  356                   if (_biharmonic._degenerate)
 
  358                       ri       -= (dphi[i][qp]*grad_u_qp)*M_prime_qp*(_biharmonic._kappa*Laplacian_u_qp);
 
  359                       ri_old   -= (dphi[i][qp]*grad_u_old_qp)*M_prime_old_qp*(_biharmonic._kappa*Laplacian_u_old_qp);
 
  362                   if (_biharmonic._cahn_hillard)
 
  366                           ri += Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp;
 
  367                           ri_old += Laplacian_phi_qp[i]*M_old_qp*_biharmonic._theta_c*(u_old_qp*u_old_qp - 1.0)*u_old_qp;
 
  368                           if (_biharmonic._degenerate)
 
  370                               ri     += (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp;
 
  371                               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;
 
  377                           ri -= Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*u_qp;
 
  378                           ri_old -= Laplacian_phi_qp[i]*M_old_qp*_biharmonic._theta_c*u_old_qp;
 
  379                           if (_biharmonic._degenerate)
 
  381                               ri     -= (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*u_qp;
 
  382                               ri_old -= (dphi[i][qp]*grad_u_old_qp)*M_prime_old_qp*_biharmonic._theta_c*u_old_qp;
 
  388                           switch(_biharmonic._log_truncation)
 
  399                   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));
 
  405                   Number M_prime_prime_qp = 0.0;
 
  406                   if (_biharmonic._degenerate) M_prime_prime_qp = -2.0;
 
  407                   for (
unsigned int j=0; j<
n_dofs; j++)
 
  410                       ri_j -= Laplacian_phi_qp[i]*M_qp*_biharmonic._kappa*Laplacian_phi_qp[j];
 
  411                       if (_biharmonic._degenerate)
 
  414                             Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._kappa*Laplacian_u_qp               +
 
  415                             (dphi[i][qp]*dphi[j][qp])*M_prime_qp*(_biharmonic._kappa*Laplacian_u_qp)                  +
 
  416                             (dphi[i][qp]*grad_u_qp)*(M_prime_prime_qp*phi[j][qp])*(_biharmonic._kappa*Laplacian_u_qp) +
 
  417                             (dphi[i][qp]*grad_u_qp)*(M_prime_qp)*(_biharmonic._kappa*Laplacian_phi_qp[j]);
 
  420                       if (_biharmonic._cahn_hillard)
 
  425                                 Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp +
 
  426                                 Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*(3.0*u_qp*u_qp - 1.0)*phi[j][qp]        +
 
  427                                 (dphi[i][qp]*dphi[j][qp])*M_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp      +
 
  428                                 (dphi[i][qp]*grad_u_qp)*M_prime_prime_qp*_biharmonic._theta_c*(u_qp*u_qp - 1.0)*u_qp  +
 
  429                                 (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*(3.0*u_qp*u_qp - 1.0)*phi[j][qp];
 
  435                                 Laplacian_phi_qp[i]*M_prime_qp*phi[j][qp]*_biharmonic._theta_c*u_qp                   +
 
  436                                 Laplacian_phi_qp[i]*M_qp*_biharmonic._theta_c*phi[j][qp]                              +
 
  437                                 (dphi[i][qp]*dphi[j][qp])*M_prime_qp*_biharmonic._theta_c*u_qp                        +
 
  438                                 (dphi[i][qp]*grad_u_qp)*M_prime_prime_qp*_biharmonic._theta_c*u_qp                    +
 
  439                                 (dphi[i][qp]*grad_u_qp)*M_prime_qp*_biharmonic._theta_c*phi[j][qp];
 
  444                               switch(_biharmonic._log_truncation)
 
  455                       Je(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] - 0.5*_biharmonic._dt*(2.0-_biharmonic._cnWeight)*ri_j);
 
  481 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 
  498   PerfLog perf_log (
"Biharmonic bounds", 
false);
 
  504   const DofMap & dof_map = get_dof_map();
 
  514   std::unique_ptr<FEBase> fe (FEBase::build(_biharmonic._dim, fe_type));
 
  522   std::vector<dof_id_type> dof_indices;
 
  524   for (
const auto & elem : _biharmonic._mesh.active_local_element_ptr_range())
 
  527       const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  534       const unsigned int n_dofs =
 
  535         cast_int<unsigned int>(dof_indices.size());
 
  542       std::vector<Point> nodes;
 
  543       fe->get_refspace_nodes(elem->type(), nodes);
 
  546       fe->reinit(elem, &nodes);
 
  559       for (
unsigned int i=0; i<
n_dofs; ++i)
 
  562           Real infinity = 1.0e20;
 
  563           Real bound = infinity;
 
  564           for (
unsigned int j = 0; j < nodes.size(); ++j)
 
  582       XL.
insert(XLe, dof_indices);
 
  583       XU.
insert(XUe, dof_indices);