Hi,
  I had a look at introduction_ex4. I found with a DirichletBoundary
object, the diagonal entries in the matrix for boundary nodes were modified
to 1 or 2. I want to add this method to adaptivity_ex2.
  It looks quite easy. I just created a DirichletBoundary object and
attached it to the system. I commented out some code with #if 0 and ran it
with *mpirun -n 1 ./.libs/example-dbg -n_timesteps 25 -n_refinements 0
-output_freq 10 -init_timestep 0*
  However, I found the dirichlet_bc object did not take effect. After some
debugging, I found libmesh somehow did not recognize boundary edges and
nodes.  I don't know what is wrong.
  I attach my modified adaptivity_ex2.C and the diff to see if someone can
help on that.
  Thank you!


--Junchao Zhang
diff --git a/examples/adaptivity/adaptivity_ex2/adaptivity_ex2.C 
b/examples/adaptivity/adaptivity_ex2/adaptivity_ex2.C
index 741b29fb..06e84767 100644
--- a/examples/adaptivity/adaptivity_ex2/adaptivity_ex2.C
+++ b/examples/adaptivity/adaptivity_ex2/adaptivity_ex2.C
@@ -51,6 +51,8 @@
 #include "libmesh/numeric_vector.h"
 #include "libmesh/dense_matrix.h"
 #include "libmesh/dense_vector.h"
+#include "libmesh/dirichlet_boundaries.h"
+#include "libmesh/analytic_function.h"
 
 #include "libmesh/getpot.h"
 
@@ -105,7 +107,13 @@ Number exact_value (const Point & p,
   return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
 }
 
-
+// Define a wrapper for exact_solution that will be needed below
+void exact_solution_wrapper (DenseVector<Number>& output,
+                             const Point& p,
+                             const Real t)
+{
+  output(0) = exact_solution(p(0), p(1), t);
+}
 
 // Begin the main program.  Note that the first
 // statement in the program throws an error if
@@ -196,9 +204,11 @@ int main (int argc, char ** argv)
   EquationSystems equation_systems (mesh);
   MeshRefinement mesh_refinement (mesh);
 
+#if 0
   // First we process the case where we do not read in the solution
   if (!read_solution)
     {
+#endif
       // Read the mesh from file.
       mesh.read ("mesh.xda");
 
@@ -222,15 +232,42 @@ int main (int argc, char ** argv)
 
       // Adds the variable "u" to "Convection-Diffusion".  "u"
       // will be approximated using first-order approximation.
-      system.add_variable ("u", FIRST);
+      unsigned int u_var = system.add_variable ("u", FIRST);
 
       // Give the system a pointer to the matrix assembly
       // and initialization functions.
       system.attach_assemble_function (assemble_cd);
       system.attach_init_function (init_cd);
 
+      int dim = 2;
+      std::set<boundary_id_type> boundary_ids;
+      boundary_ids.insert(0);
+      boundary_ids.insert(1);
+      if(dim>=2)
+        {
+          boundary_ids.insert(2);
+          boundary_ids.insert(3);
+        }
+      if(dim==3)
+        {
+          boundary_ids.insert(4);
+          boundary_ids.insert(5);
+        }
+
+      std::vector<unsigned int> variables(1);
+      variables[0] = u_var;
+
+      AnalyticFunction<> exact_solution_object(exact_solution_wrapper);
+
+      DirichletBoundary dirichlet_bc(boundary_ids,
+                                     variables,
+                                     &exact_solution_object);
+
+      system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
+
       // Initialize the data structures for the equation system.
       equation_systems.init ();
+#if 0
     }
   // Otherwise we read in the solution and mesh
   else
@@ -262,6 +299,7 @@ int main (int argc, char ** argv)
 
       libMesh::out << "Initial H1 norm = " << H1norm << std::endl << std::endl;
     }
+#endif
 
   // Prints information about the system to the screen.
   equation_systems.print_info();
@@ -296,8 +334,10 @@ int main (int argc, char ** argv)
   // derived from them) contain old solutions, to use the
   // old_local_solution later we now need to specify the system
   // type when we ask for it.
+#if 0
   TransientLinearImplicitSystem & system =
     
equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
+#endif
 
   const Real dt = 0.025;
   system.time = init_timestep*dt;
@@ -645,6 +685,7 @@ void assemble_cd (EquationSystems & es,
             }
         }
 
+#if 0
       // At this point the interior element integration has
       // been completed.  However, we have not yet addressed
       // boundary conditions.  For this example we will only
@@ -683,6 +724,7 @@ void assemble_cd (EquationSystems & es,
                 }
             }
       }
+#endif
 
 
       // We have now built the element matrix and RHS vector in terms
------------------------------------------------------------------------------
Site24x7 APM Insight: Get Deep Visibility into Application Performance
APM + Mobile APM + RUM: Monitor 3 App instances at just $35/Month
Monitor end-to-end web transactions and take corrective actions now
Troubleshoot faster and improve end-user experience. Signup Now!
http://pubads.g.doubleclick.net/gampad/clk?id=267308311&iu=/4140
_______________________________________________
Libmesh-users mailing list
Libmesh-users@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/libmesh-users

Reply via email to