]> https://gitweb.dealii.org/ - code-gallery.git/commitdiff
Some fixes for functions of make_grid, make_constraints, and print_vertical_tip_displ...
authorlianxiali <lianxiali@gmail.com>
Mon, 26 Jul 2021 03:31:26 +0000 (21:31 -0600)
committerlianxiali <lianxiali@gmail.com>
Mon, 26 Jul 2021 03:31:26 +0000 (21:31 -0600)
Quasi_static_Finite_strain_Compressible_Elasticity/cook_membrane.cc

index 9b19f0c68f25b8ee9a2db9726206568f5e65963c..f95f8390235c3942ce9f9ca3eecebeed8147b681 100644 (file)
@@ -1101,7 +1101,7 @@ namespace Cook_Membrane
               cell->face(face)->set_boundary_id(1); // -X faces
             else if (std::abs(cell->face(face)->center()[0] - 48.0) < tol_boundary)
               cell->face(face)->set_boundary_id(11); // +X faces
-            else if (std::abs(std::abs(cell->face(face)->center()[0]) - 0.5) < tol_boundary)
+            else if (std::abs(std::abs(cell->face(face)->center()[2]) - 0.5) < tol_boundary)
               cell->face(face)->set_boundary_id(2); // +Z and -Z faces
           }
 
@@ -1369,7 +1369,9 @@ namespace Cook_Membrane
       std::cout << "_";
     std::cout << std::endl;
 
-    Point<dim> soln_pt (48.0*parameters.scale,60.0*parameters.scale);
+    Point<dim> soln_pt{};
+    soln_pt[0] = 48.0*parameters.scale;
+    soln_pt[1] = 60.0*parameters.scale;
     if (dim == 3)
       soln_pt[2] = 0.5*parameters.scale;
     double vertical_tip_displacement = 0.0;
@@ -2062,12 +2064,12 @@ namespace Cook_Membrane
     const bool apply_dirichlet_bc = (it_nr == 0);
 
     // The boundary conditions for the indentation problem are as follows: On
-    // the -x, -y and -z faces (ID's 0,2,4) we set up a symmetry condition to
-    // allow only planar movement while the +x and +y faces (ID's 1,3) are
-    // traction free. In this contrived problem, part of the +z face (ID 5) is
-    // set to have no motion in the x- and y-component. Finally, as described
-    // earlier, the other part of the +z face has an the applied pressure but
-    // is also constrained in the x- and y-directions.
+    // the -x face (ID = 1), we set up a zero-displacement condition, -y and +y traction 
+    // free boundary condition (don't need to take care); -z and +z faces (ID = 2) are 
+    // not allowed to move along z axis so that it is a plane strain problem. 
+    // Finally, as described earlier, +x face (ID = 11) has an the applied 
+    // distributed shear force (converted by total force per unit area) which 
+    // needs to be taken care as an inhomogeneous Newmann boundary condition.
     //
     // In the following, we will have to tell the function interpolation
     // boundary values which components of the solution vector should be
@@ -2088,12 +2090,6 @@ namespace Cook_Membrane
                                                  ZeroFunction<dim>(n_components),
                                                  constraints,
                                                  fe.component_mask(u_fe));
-      else
-        VectorTools::interpolate_boundary_values(dof_handler_ref,
-                                                 boundary_id,
-                                                 ZeroFunction<dim>(n_components),
-                                                 constraints,
-                                                 fe.component_mask(u_fe));
     }
 
     // Zero Z-displacement through thickness direction
@@ -2109,12 +2105,6 @@ namespace Cook_Membrane
                                                    ZeroFunction<dim>(n_components),
                                                    constraints,
                                                    fe.component_mask(z_displacement));
-        else
-          VectorTools::interpolate_boundary_values(dof_handler_ref,
-                                                   boundary_id,
-                                                   ZeroFunction<dim>(n_components),
-                                                   constraints,
-                                                   fe.component_mask(z_displacement));
       }
 
     constraints.close();

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.