84 * #include <deal.II/base/quadrature_lib.h>
85 * #include <deal.II/base/function.h>
86 * #include <deal.II/base/logstream.h>
87 * #include <deal.II/base/utilities.h>
88 * #include <deal.II/base/convergence_table.h>
89 * #include <deal.II/base/smartpointer.h>
90 * #include <deal.II/base/parameter_handler.h>
91 * #include <deal.II/base/timer.h>
93 * #include <deal.II/lac/vector.h>
94 * #include <deal.II/lac/full_matrix.h>
95 * #include <deal.II/lac/sparse_matrix.h>
96 * #include <deal.II/lac/solver_cg.h>
97 * #include <deal.II/lac/precondition.h>
98 * #include <deal.II/lac/affine_constraints.h>
99 * #include <deal.II/lac/dynamic_sparsity_pattern.h>
101 * #include <deal.II/grid/
tria.h>
102 * #include <deal.II/grid/grid_generator.h>
103 * #include <deal.II/grid/tria_accessor.h>
104 * #include <deal.II/grid/tria_iterator.h>
105 * #include <deal.II/grid/manifold_lib.h>
106 * #include <deal.II/grid/grid_refinement.h>
107 * #include <deal.II/grid/grid_in.h>
109 * #include <deal.II/dofs/dof_handler.h>
110 * #include <deal.II/dofs/dof_accessor.h>
111 * #include <deal.II/dofs/dof_tools.h>
112 * #include <deal.II/dofs/dof_renumbering.h>
114 * #include <deal.II/fe/fe_values.h>
115 * #include <deal.II/fe/fe_q.h>
117 * #include <deal.II/numerics/vector_tools.h>
118 * #include <deal.II/numerics/matrix_tools.h>
119 * #include <deal.II/numerics/data_out.h>
120 * #include <deal.II/numerics/error_estimator.h>
121 * #include <deal.II/numerics/solution_transfer.h>
123 * #include <typeinfo>
125 * #include <iostream>
127 * #include <deal.II/numerics/solution_transfer.h>
131 * Open a
namespace for this program and import everything from the
132 *
dealii namespace into it.
141 * ********************************************************
148 *
void read_parameters(
const std::string);
150 *
void declare_parameters();
162 *
void ParameterReader::declare_parameters()
165 * prm.enter_subsection (
"Global Parameters");
168 *
"Penalization parameter");
170 *
"Whether the exact solution is known");
172 * prm.leave_subsection ();
174 * prm.enter_subsection (
"Mesh & Refinement Parameters");
177 *
"Number identifying the domain in which we solve the problem");
179 *
"Number of global mesh refinement steps applied to initial coarse grid");
181 *
"Number of global adaptive mesh refinements");
183 *
"refinement threshold");
185 *
"coarsening threshold");
187 * prm.leave_subsection ();
190 * prm.enter_subsection (
"Algorithm Parameters");
193 *
"0: Preconditioned descent, 1: Newton Method");
199 *
"Maximum Number of CG iterations");
201 *
"Tolerance for CG iterations");
203 *
"Maximum Number of LS iterations");
205 *
"line search tolerance constant (c1 in Nocedal-Wright)");
207 *
"initial step length in line-search");
209 *
"Maximum Number of inner iterations");
211 *
"Threshold on norm of the derivative to declare optimality achieved");
213 *
"Threshold on norm of the derivative to declare optimality achieved in highly refined mesh");
215 *
"Number of adaptive refinement before change convergence threshold");
217 * prm.leave_subsection ();
220 *
void ParameterReader::read_parameters (
const std::string parameter_file)
222 * declare_parameters();
223 * prm.parse_input (parameter_file);
228 * ******************************************************************************************
229 * The solution of the elastoplastic torsion problem on the unit disk with rhs=4.
236 *
class Solution :
public Function<dim>
240 *
virtual double value (
const Point<dim> &pto,
const unsigned int component = 0)
const;
245 *
double Solution<dim>::value (
const Point<dim> &pto,
const unsigned int)
const
247 *
double r=
sqrt(pto.square());
259 *
double r=
sqrt(pto.square());
271 * ******************************************************************************************
283 * ComputeMultiplier (
double pe);
286 *
void compute_derived_quantities_scalar (
287 *
const std::vector< double > &,
295 *
virtual std::vector<std::string>
get_names ()
const;
298 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
305 * ComputeMultiplier<dim>::ComputeMultiplier (
double pe): p(pe)
310 *
void ComputeMultiplier<dim>::compute_derived_quantities_scalar(
311 *
const std::vector< double > & ,
318 *
const unsigned int n_quadrature_points = duh.size();
320 *
for (
unsigned int q=0; q<n_quadrature_points; ++q)
322 *
long double sqrGrad=duh[q]* duh[q];
323 *
long double exponent=(p-2.0)/2*
std::log(sqrGrad);
324 * computed_quantities[q](0) =
std::sqrt(sqrGrad);
325 * computed_quantities[q](1)=
std::exp(exponent);
334 * std::vector<std::string>
335 * ComputeMultiplier<dim>::get_names() const
337 * std::vector<std::string> solution_names;
338 * solution_names.push_back (
"Gradient norm");
339 * solution_names.push_back (
"Lagrange multiplier");
340 *
return solution_names;
346 * ComputeMultiplier<dim>::get_needed_update_flags () const
354 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
355 * ComputeMultiplier<dim>:: get_data_component_interpretation () const
357 * std::vector<DataComponentInterpretation::DataComponentInterpretation>
367 * Lagrange multiplier
371 *
return interpretation;
380 * ***************************************************************************************
384 *
class ElastoplasticTorsion
388 * ~ElastoplasticTorsion ();
392 *
void setup_system (
const bool initial_step);
393 *
void assemble_system ();
394 *
bool solve (
const int inner_it);
396 *
void refine_mesh ();
397 *
void set_boundary_values ();
398 *
double phi (
const double alpha)
const;
399 *
bool checkWolfe(
double &alpha,
double &phi_alpha)
const;
400 *
bool determine_step_length (
const int inner_it);
401 *
void print_it_message (
const int counter,
bool ks);
402 *
void output_results (
unsigned int refinement)
const;
403 *
void format_convergence_tables();
404 *
void process_solution (
const unsigned int cycle);
405 *
void process_multiplier (
const unsigned int cycle,
const int iter,
double time);
406 *
double dual_error ()
const;
407 *
double dual_infty_error ()
const;
408 *
double W (
double Du2)
const;
409 *
double Wp (
double Du2)
const;
410 *
double G (
double Du2)
const;
429 *
double step_length,phi_zero,phi_alpha,phip,phip_zero;
430 *
double old_step,old_phi_zero,old_phip;
433 *
double Linfty_error;
434 *
double dual_L1_error;
435 *
double dual_L_infty_error;
438 *
double line_search_tolerence;
439 *
unsigned int dir_id;
440 * std::string elements;
441 * std::string Method;
455 *
class BoundaryValues :
public Function<dim>
458 * BoundaryValues () :
Function<dim>() {}
461 *
const unsigned int component = 0)
const;
466 *
double BoundaryValues<dim>::value (
const Point<dim> &,
467 *
const unsigned int )
const
471 * could be anything
else (theory works provided |Dg|_infty < 1/2)
505 *
class RightHandSide :
public Function<dim>
508 * RightHandSide () :
Function<dim>() {}
510 *
const unsigned int component = 0)
const;
514 *
double RightHandSide<dim>::value (
const Point<dim> &,
515 *
const unsigned int )
const
519 *
set to
constant = 4,
for which
explicit solution to compare exists
523 *
double return_value = 4.0;
524 *
return return_value;
532 * The ElastoplasticTorsion
class implementation
536 * Constructor of the
class
540 * ElastoplasticTorsion<dim>::ElastoplasticTorsion (
ParameterHandler ¶m):
546 * dual_L1_error(1.0),
547 * dual_L_infty_error(1.0),
550 * prm.enter_subsection (
"Global Parameters");
551 * p=prm.get_double(
"p");
552 * prm.leave_subsection ();
553 * prm.enter_subsection (
"Algorithm Parameters");
554 * line_search_tolerence=prm.get_double(
"line_search_tolerence");
555 * dir_id=prm.get_integer(
"Descent_direction");
556 * prm.leave_subsection ();
559 *
else elements=
"P2";
570 * ElastoplasticTorsion<dim>::~ElastoplasticTorsion ()
572 * dof_handler.clear ();
578 * print iteration message
585 *
void ElastoplasticTorsion<dim>::print_it_message (
const int counter,
bool ks)
589 * process_solution (counter);
590 * std::cout <<
"iteration="<< counter+1 <<
" J(u_h)= "<< phi_zero <<
", H1 error: "
591 * << H1_error <<
", W0-1,infty error: "<< Linfty_error<<
" J'(u_h)(w)= "<< phip
592 * <<
", |J'(u_h)|= "<< system_rhs.l2_norm()<<std::endl;
596 * std::cout <<
"iteration= " << counter+1 <<
" J(u_h)= "
597 * << phi_alpha <<
" J'(u_h)= "<< phip<<std::endl;
624 *
void ElastoplasticTorsion<dim>::format_convergence_tables()
626 * convergence_table.set_precision(
"L2", 3);
627 * convergence_table.set_precision(
"H1", 3);
628 * convergence_table.set_precision(
"Linfty", 3);
629 * convergence_table.set_precision(
"function value", 3);
630 * convergence_table.set_precision(
"derivative", 3);
631 * dual_convergence_table.set_precision(
"dual_L1", 3);
632 * dual_convergence_table.set_precision(
"dual_Linfty", 3);
633 * dual_convergence_table.set_precision(
"L2", 3);
634 * dual_convergence_table.set_precision(
"H1", 3);
635 * dual_convergence_table.set_precision(
"Linfty", 3);
636 * convergence_table.set_scientific(
"L2",
true);
637 * convergence_table.set_scientific(
"H1",
true);
638 * convergence_table.set_scientific(
"Linfty",
true);
639 * convergence_table.set_scientific(
"function value",
true);
640 * convergence_table.set_scientific(
"derivative",
true);
641 * dual_convergence_table.set_scientific(
"dual_L1",
true);
642 * dual_convergence_table.set_scientific(
"dual_Linfty",
true);
643 * dual_convergence_table.set_scientific(
"L2",
true);
644 * dual_convergence_table.set_scientific(
"H1",
true);
645 * dual_convergence_table.set_scientific(
"Linfty",
true);
652 * fill-in entry
for the solution
656 *
void ElastoplasticTorsion<dim>::process_solution (
const unsigned int it)
662 * compute
L2 error (save to difference_per_cell)
667 * L2_error = difference_per_cell.l2_norm();
671 * compute H1 error (save to difference_per_cell)
676 * H1_error = difference_per_cell.l2_norm();
680 * compute W1infty error (save to difference_per_cell)
687 * Linfty_error = difference_per_cell.linfty_norm();
690 * convergence_table.add_value(
"cycle", it);
691 * convergence_table.add_value(
"p", p);
692 * convergence_table.add_value(
"L2", L2_error);
693 * convergence_table.add_value(
"H1", H1_error);
694 * convergence_table.add_value(
"Linfty", Linfty_error);
695 * convergence_table.add_value(
"function value", phi_alpha);
696 * convergence_table.add_value(
"derivative", phip);
703 * fill-in entry
for the multiplier
707 *
void ElastoplasticTorsion<dim>::process_multiplier (
const unsigned int cycle,
const int iter,
double time)
710 *
const unsigned int n_dofs=dof_handler.n_dofs();
711 * dual_L1_error=dual_error();
712 * dual_L_infty_error=dual_infty_error();
715 * dual_convergence_table.add_value(
"cycle", cycle);
716 * dual_convergence_table.add_value(
"p", p);
717 * dual_convergence_table.add_value(
"iteration_number", iter);
718 * dual_convergence_table.add_value(
"cpu_time", time);
719 * dual_convergence_table.add_value(
"cells", n_active_cells);
720 * dual_convergence_table.add_value(
"dofs", n_dofs);
721 * dual_convergence_table.add_value(
"L2", L2_error);
722 * dual_convergence_table.add_value(
"H1", H1_error);
723 * dual_convergence_table.add_value(
"Linfty", Linfty_error);
724 * dual_convergence_table.add_value(
"dual_L1", dual_L1_error);
725 * dual_convergence_table.add_value(
"dual_Linfty", dual_L_infty_error);
735 * ElastoplasticTorsion::setup_system
736 * unchanged from @ref step_15
"step-15"
743 *
void ElastoplasticTorsion<dim>::setup_system (
const bool initial_step)
747 * dof_handler.distribute_dofs (fe);
748 * present_solution.reinit (dof_handler.n_dofs());
749 * grad_norm.reinit (dof_handler.n_dofs());
750 *
lambda.reinit (dof_handler.n_dofs());
752 * hanging_node_constraints.clear ();
754 * hanging_node_constraints);
755 * hanging_node_constraints.close ();
761 * The remaining parts of the function
767 * newton_update.reinit (dof_handler.n_dofs());
768 * system_rhs.reinit (dof_handler.n_dofs());
771 * hanging_node_constraints.condense (c_sparsity);
772 * sparsity_pattern.copy_from(c_sparsity);
773 * system_matrix.reinit (sparsity_pattern);
785 *
double ElastoplasticTorsion<dim>::W (
double Du2)
const
791 *
double ElastoplasticTorsion<dim>::Wp (
double )
const
797 *
double ElastoplasticTorsion<dim>::G (
double )
const
804 *
void ElastoplasticTorsion<dim>::assemble_system ()
809 *
const RightHandSide<dim> right_hand_side;
819 *
const unsigned int dofs_per_cell = fe.dofs_per_cell;
820 *
const unsigned int n_q_points = quadrature_formula.size();
825 * std::vector<Tensor<1, dim> > old_solution_gradients(n_q_points);
826 * std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
830 * cell = dof_handler.begin_active(),
831 * endc = dof_handler.end();
832 *
for (; cell!=endc; ++cell)
837 * fe_values.reinit (cell);
838 * fe_values.get_function_gradients(present_solution,
839 * old_solution_gradients);
841 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
843 *
long double coeff=0.0;
844 *
long double a=old_solution_gradients[q_point] * old_solution_gradients[q_point];
845 *
long double exponent=(p-2.0)/2*
std::log(a);
847 *
for (
unsigned int i=0; i<dofs_per_cell; ++i)
849 *
for (
unsigned int j=0; j<dofs_per_cell; ++j)
853 *
cell_matrix(i, j) += fe_values.shape_grad(i, q_point) * fe_values.shape_grad(j, q_point)
854 * * (G(a)+(p-1.0)*coeff) * fe_values.JxW(q_point);
858 *
cell_matrix(i, j) += fe_values.shape_grad(i, q_point) * fe_values.shape_grad(j, q_point)
860 * * fe_values.JxW(q_point);
864 * cell_rhs(i) -= ( fe_values.shape_grad(i, q_point)
865 * * old_solution_gradients[q_point]
867 * -right_hand_side.value(fe_values.quadrature_point(q_point))
868 * *fe_values.shape_value(i, q_point)
870 * * fe_values.JxW(q_point);
874 * cell->get_dof_indices (local_dof_indices);
875 *
for (
unsigned int i=0; i<dofs_per_cell; ++i)
877 *
for (
unsigned int j=0; j<dofs_per_cell; ++j)
878 * system_matrix.add (local_dof_indices[i],
879 * local_dof_indices[j],
882 * system_rhs(local_dof_indices[i]) += cell_rhs(i);
886 * hanging_node_constraints.condense (system_matrix);
887 * hanging_node_constraints.condense (system_rhs);
889 * std::map<types::global_dof_index,double> boundary_values;
906 * unchanged from @ref step_15
"step-15"
913 *
void ElastoplasticTorsion<dim>::refine_mesh ()
915 *
using FunctionMap = std::map<types::boundary_id, const Function<dim> *>;
922 * estimated_error_per_cell);
924 * prm.enter_subsection (
"Mesh & Refinement Parameters");
925 *
const double top_fraction=prm.get_double(
"top_fraction_of_cells");
926 *
const double bottom_fraction=prm.get_double(
"bottom_fraction_of_cells");
927 * prm.leave_subsection ();
929 * estimated_error_per_cell,
930 * top_fraction, bottom_fraction);
934 * solution_transfer.prepare_for_coarsening_and_refinement(present_solution);
936 * dof_handler.distribute_dofs(fe);
938 * solution_transfer.interpolate(present_solution, tmp);
939 * present_solution = tmp;
940 * set_boundary_values ();
941 * hanging_node_constraints.clear();
944 * hanging_node_constraints);
945 * hanging_node_constraints.close();
946 * hanging_node_constraints.distribute (present_solution);
947 * setup_system (
false);
954 * Dump the
norm of the
gradient and the lagrange multiplier in
vtu format
for visualization
958 *
void ElastoplasticTorsion<dim>::output_results (
unsigned int counter)
const
962 * multiplier
object contains both |Du| and
lambda.
965 * ComputeMultiplier<dim> multiplier(p);
972 * std::ostringstream p_str;
973 * p_str << p<<
"-cycle-"<<counter;
974 * std::string str = p_str.str();
975 *
const std::string filename =
"solution-" + str+
".vtu";
976 * std::ofstream output (filename.c_str());
983 * unchanged from @ref step_15
"step-15"
987 *
void ElastoplasticTorsion<dim>::set_boundary_values ()
989 * std::map<types::global_dof_index, double> boundary_values;
992 * BoundaryValues<dim>(),
994 *
for (std::map<types::global_dof_index, double>::const_iterator
995 * bp = boundary_values.begin();
996 * bp != boundary_values.end(); ++bp)
997 * present_solution(bp->first) = bp->second;
1004 * COMPUTE @f$\phi(\alpha)=J_p(u_h+\alpha w)@f$
1007 *
template <
int dim>
1008 *
double ElastoplasticTorsion<dim>::phi (
const double alpha)
const
1011 *
const RightHandSide<dim> right_hand_side;
1013 * evaluation_point = present_solution;
1014 * evaluation_point.add (alpha, newton_update);
1023 *
const unsigned int dofs_per_cell = fe.dofs_per_cell;
1024 *
const unsigned int n_q_points = quadrature_formula.size();
1027 * std::vector<Tensor<1, dim> >
gradients(n_q_points);
1028 * std::vector<double>
values(n_q_points);
1031 * std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1035 * endc = dof_handler.end();
1036 *
for (; cell!=endc; ++cell)
1039 * fe_values.reinit (cell);
1040 * fe_values.get_function_gradients (evaluation_point, gradients);
1041 * fe_values.get_function_values (evaluation_point, values);
1044 *
for (
unsigned int q_point=0; q_point<n_q_points; ++q_point)
1055 * obj+= 1/2 W(|Du|^2)+1/p |Du|^p -fu (see (1))
1059 * (0.5*W(Du2)+penalty/p)- right_hand_side.value(fe_values.quadrature_point(q_point))*values[q_point]
1060 * ) * fe_values.JxW(q_point);
1072 * Compute L^1 error
norm of Lagrange Multiplier
1073 * with respect to exact solution (cf. Alvarez & Flores, 2015)
1079 *
template <
int dim>
1080 *
double ElastoplasticTorsion<dim>::dual_error () const
1090 *
const unsigned int dofs_per_cell = fe.dofs_per_cell;
1091 *
const unsigned int n_q_points = quadrature_formula.size();
1094 * std::vector<Tensor<1, dim> >
gradients(n_q_points);
1096 * std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1100 * endc = dof_handler.end();
1101 *
for (; cell!=endc; ++cell)
1104 * fe_values.reinit (cell);
1105 * fe_values.get_function_gradients (present_solution, gradients);
1107 *
for (
unsigned int q_point=0; q_point<n_q_points; ++q_point)
1110 *
if (coeff<1.0e-15)
1113 * coeff=
std::pow(coeff,(p-2.0)/2.0);
1115 *
double r=
std::sqrt(fe_values.quadrature_point(q_point).square());
1120 * obj+=(
std::abs(coeff-exact) ) * fe_values.JxW(q_point);
1131 * Compute L^infinity error
norm of Lagrange Multiplier
1132 * with respect to exact solution (cf. Alvarez & Flores, 2015)
1138 *
template <
int dim>
1139 *
double ElastoplasticTorsion<dim>::dual_infty_error () const
1149 *
const unsigned int dofs_per_cell = fe.dofs_per_cell;
1150 *
const unsigned int n_q_points = quadrature_formula.size();
1153 * std::vector<Tensor<1, dim> >
gradients(n_q_points);
1155 * std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1159 * endc = dof_handler.end();
1160 *
for (; cell!=endc; ++cell)
1163 * fe_values.reinit (cell);
1164 * fe_values.get_function_gradients (present_solution, gradients);
1166 *
for (
unsigned int q_point=0; q_point<n_q_points; ++q_point)
1169 *
double r=
std::sqrt(fe_values.quadrature_point(q_point).square());
1175 * compute |Du|^(p-2) as
exp(p-2/2*
log(Du^2))
1178 *
long double exponent=(p-2.0)/2*
std::log(sqdGrad);
1179 *
long double coeff=
std::exp(exponent);
1193 *
check whether putative step-length satisfies sufficient decrease conditions
1196 *
template <
int dim>
1197 *
bool ElastoplasticTorsion<dim>::checkWolfe(
double &alpha,
double &phi_alpha)
const
1199 *
if (phi_alpha< phi_zero+line_search_tolerence*phip*alpha )
1209 * Find a step-length satisfying sufficient decrease condition by line-search
1210 * uses quadratic interpolation
1216 *
template <
int dim>
1217 *
bool ElastoplasticTorsion<dim>::determine_step_length(
const int inner_it)
1219 *
unsigned int it=0;
1221 *
double alpha,nalpha;
1222 * prm.enter_subsection (
"Algorithm Parameters");
1223 *
const unsigned int max_LS_it=prm.get_integer(
"max_LS_it");
1224 *
double init_SL=prm.get_double(
"init_step_length");
1225 * prm.leave_subsection ();
1230 * alpha=
std::min(1.45*old_step*old_phip/phip,1.0);
1232 * phi_alpha=phi(alpha);
1233 * std::cerr <<
"Step length=" << alpha <<
", Value= " << phi_alpha;
1236 *
check if step-size satisfies sufficient decrease condition
1239 * done=checkWolfe(alpha,phi_alpha);
1241 * std::cerr <<
" accepted" << std::endl;
1243 * std::cerr <<
" rejected" ;
1245 *
while ((!done) & (it<max_LS_it))
1249 *
new try obtained by quadratic interpolation
1252 * nalpha=-(phip*alpha*alpha)/(2*(phi_alpha-phi_zero-phip*alpha));
1254 *
if (nalpha<1e-3*alpha ||
std::abs(nalpha-alpha)/alpha<1
e-8)
1256 *
else if ( phi_alpha-phi_zero>1e3*
std::abs(phi_zero) )
1259 * phi_alpha=phi(alpha);
1260 * done=checkWolfe(alpha,phi_alpha);
1262 * std::cerr <<
", finished with steplength= "<< alpha<<
", fcn value= "<< phi_alpha<<std::endl;
1267 * std::cerr <<
", max. no. of iterations reached wiht steplength= "<< alpha
1268 * <<
", fcn value= "<< phi_alpha<<std::endl;
1273 * step_length=alpha;
1282 * ElastoplasticTorsion::init_mesh()
1288 *
template <
int dim>
1289 *
void ElastoplasticTorsion<dim>::init_mesh ()
1296 * prm.enter_subsection (
"Mesh & Refinement Parameters");
1297 *
const int domain_id=prm.get_integer(
"Code for the domain");
1298 *
const int init_ref=prm.get_integer(
"No of initial refinements");
1299 * prm.leave_subsection ();
1306 * For the unit disk around the origin
1313 *
else if (domain_id==1)
1317 * For the unit square
1322 *
else if (domain_id==2)
1362 * ElastoplasticTorsion::solve(inner_it)
1363 * Performs one inner iteration
1369 *
template <
int dim>
1370 *
bool ElastoplasticTorsion<dim>::solve (
const int inner_it)
1372 * prm.enter_subsection (
"Algorithm Parameters");
1373 *
const unsigned int max_CG_it=prm.get_integer(
"Max_CG_it");
1374 *
const double CG_tol=prm.get_double(
"CG_tol");
1375 * prm.leave_subsection ();
1381 * preconditioner.
initialize(system_matrix,0.25);
1383 * solver.solve (system_matrix, newton_update, system_rhs,
1385 * hanging_node_constraints.distribute (newton_update);
1389 * Recall that phi(alpha)=J(u+alpha w)
1392 * old_step=step_length;
1393 * old_phi_zero=phi_zero;
1396 * phip=-1.0*(newton_update*system_rhs);
1402 * std::cout <<
"Not a descent direction!" <<std::endl;
1403 * present_solution.add (-1.0*step_length, newton_update);
1404 * step_length=step_length/2;
1410 *
if (determine_step_length(inner_it))
1414 * update u_{n+1}=u_n+alpha w_n
1417 * present_solution.add (step_length, newton_update);
1420 *
else return false;
1429 * ElastoplasticTorsion::run
1432 *
template <
int dim>
1433 *
void ElastoplasticTorsion<dim>::run ()
1441 * prm.enter_subsection (
"Mesh & Refinement Parameters");
1442 *
const int adapt_ref=prm.get_integer(
"No of adaptive refinements");
1443 * prm.leave_subsection ();
1444 * prm.enter_subsection (
"Algorithm Parameters");
1445 *
const int max_inner=prm.get_integer(
"Max_inner");
1446 *
const double eps=prm.get_double(
"eps");
1447 *
const double hi_eps=prm.get_double(
"hi_eps");
1448 *
const int hi_th=prm.get_integer(
"hi_th");
1449 *
const double init_p=prm.get_double(
"init_p");
1450 *
const double delta_p=prm.get_double(
"delta_p");
1451 * prm.leave_subsection ();
1452 * prm.enter_subsection (
"Global Parameters");
1453 *
bool known_solution=prm.get_bool(
"known_solution");
1454 *
double actual_p=prm.get_double(
"p");
1455 * prm.leave_subsection ();
1469 * initalize mesh
for the selected domain
1479 * setup_system (
true);
1480 * set_boundary_values ();
1489 *
int refinement = 0;
1498 *
bool well_solved=
true;
1502 *
while (p<actual_p)
1504 * std::cout <<
"--Preparing initial condition with p="<<p<<
" iter.= " << global_it<<
" .-- "<< std::endl;
1506 *
for (
int inner_iteration=0; inner_iteration<max_inner; ++inner_iteration,++global_it)
1508 * assemble_system ();
1509 * well_solved=solve (inner_iteration);
1510 * print_it_message (global_it, known_solution);
1512 * ((system_rhs.l2_norm()/
std::sqrt(system_rhs.size()) <1e-4) & (cycle<1)) |
1513 * ((system_rhs.l2_norm()/
std::sqrt(system_rhs.size()) <1e-5) & (cycle>=1)) |
1518 * ptime=timer.cpu_time();
1520 * output_results (cycle);
1522 *
if (known_solution)
1524 * process_multiplier(cycle,global_it,ptime);
1527 * dual_convergence_table.write_tex(dual_error_table_file);
1550 * std::cout <<
"============ Solving problem with p=" <<p <<
" ==================" << std::endl;
1552 *
while ((cycle<adapt_ref) & well_solved)
1560 *
for (
int inner_iteration=0; inner_iteration<max_inner; ++inner_iteration,++global_it)
1562 * assemble_system ();
1563 * well_solved=solve (inner_iteration);
1564 * print_it_message (global_it, known_solution);
1567 * ((system_rhs.l2_norm()/
std::sqrt(system_rhs.size()) < eps) & (refinement<hi_th)) |
1568 * (( system_rhs.l2_norm()/
std::sqrt (system_rhs.size()) <hi_eps) | (!well_solved))
1574 * inner iterations finished
1577 * ptime=timer.cpu_time();
1579 * output_results (cycle);
1583 * compute and display error,
if the
explicit solution is known
1586 *
if (known_solution)
1588 * process_multiplier(cycle,global_it,ptime);
1589 * std::cout <<
"finished with H1 error: " << H1_error <<
", dual error (L1): "
1590 * << dual_L1_error <<
"dual error (L infty): "<<dual_L_infty_error <<std::endl;
1605 * std::cout <<
"******** Refined mesh " << cycle <<
" ********" << std::endl;
1611 * write convergence tables to file
1614 *
if (known_solution)
1616 * format_convergence_tables();
1617 * std::string error_filename =
"error"+Method+elements+
".tex";
1618 * std::ofstream error_table_file(error_filename.c_str());
1619 * std::string dual_error_filename =
"dual_error"+Method+elements+
".tex";
1620 * std::ofstream dual_error_table_file(dual_error_filename.c_str());
1621 * convergence_table.write_tex(error_table_file);
1622 * dual_convergence_table.write_tex(dual_error_table_file);
1638 *
using namespace dealii;
1639 *
using namespace nsp;
1643 * ParameterReader param(prm);
1644 * param.read_parameters(
"EPT.prm");
1645 * ElastoplasticTorsion<2> ElastoplasticTorsionProblem(prm);
1646 * ElastoplasticTorsionProblem .run ();
1648 *
catch (std::exception &exc)
1650 * std::cerr << std::endl << std::endl
1651 * <<
"----------------------------------------------------" << std::endl;
1652 * std::cerr <<
"Exception on processing: " << std::endl
1653 * << exc.what() << std::endl
1654 * <<
"Aborting!" << std::endl
1655 * <<
"----------------------------------------------------"
1662 * std::cerr << std::endl << std::endl
1663 * <<
"----------------------------------------------------"
1665 * std::cerr <<
"Unknown exception!" << std::endl
1666 * <<
"Aborting!" << std::endl
1667 * <<
"----------------------------------------------------"
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
virtual UpdateFlags get_needed_update_flags() const =0
virtual std::vector< std::string > get_names() const =0
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
active_cell_iterator begin_active(const unsigned int level=0) const
virtual Tensor< 1, dim, RangeNumberType > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
unsigned int depth_console(const unsigned int n)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
__global__ void set(Number *val, const Number s, const size_type N)
void write_vtu(std::ostream &out) const
typename ActiveSelector::active_cell_iterator active_cell_iterator
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void hyper_ball(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
void hyper_rectangle(Triangulation< dim, spacedim > &tria, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void merge_triangulations(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result, const double duplicated_vertex_tolerance=1.0e-12, const bool copy_manifold_ids=false)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, const ArrayView< const std::vector< double > > &velocity, double factor=1.)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void run(const Iterator &begin, const typename identity< Iterator >::type &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
bool check(const ConstraintKinds kind_in, const unsigned int dim)
unsigned int n_active_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const ::Triangulation< dim, spacedim > & tria