]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
add local integrators and list namespaces in module
authorkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 9 Nov 2010 15:00:06 +0000 (15:00 +0000)
committerkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 9 Nov 2010 15:00:06 +0000 (15:00 +0000)
git-svn-id: https://svn.dealii.org/trunk@22642 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/integrators/differential.h
deal.II/include/deal.II/integrators/elasticity.h [new file with mode: 0644]
deal.II/include/deal.II/integrators/l2.h
deal.II/include/deal.II/integrators/laplace.h
deal.II/include/deal.II/integrators/local_integrators.h
deal.II/include/deal.II/integrators/maxwell.h [new file with mode: 0644]
deal.II/include/deal.II/numerics/mesh_worker.h

index 5beec751cb39c2ae8aae2d53fbb3dfcc17666ad0..00af25a290816e2d290ba88dfc91912b3d2aa21c 100644 (file)
@@ -29,6 +29,7 @@ namespace LocalIntegrators
 /**
  * @brief Local integrators related to first order differential operators and their traces.
  *
+ * @ingroup Integrators
  * @author Guido Kanschat
  * @date 2010
  */
@@ -52,6 +53,7 @@ namespace LocalIntegrators
       const unsigned int n_dofs = fe.dofs_per_cell;
       const unsigned int t_dofs = fetest.dofs_per_cell;
       AssertDimension(fecomp, dim);
+      AssertDimension(fetest.get_fe().n_components(), 1);
       AssertDimension(M.m(), t_dofs);
       AssertDimension(M.n(), n_dofs);
        
diff --git a/deal.II/include/deal.II/integrators/elasticity.h b/deal.II/include/deal.II/integrators/elasticity.h
new file mode 100644 (file)
index 0000000..f30d965
--- /dev/null
@@ -0,0 +1,82 @@
+//---------------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//---------------------------------------------------------------------------
+#ifndef __deal2__integrators_elasticity_h
+#define __deal2__integrators_elasticity_h
+
+
+#include <base/config.h>
+#include <base/exceptions.h>
+#include <base/quadrature.h>
+#include <lac/full_matrix.h>
+#include <fe/mapping.h>
+#include <fe/fe_values.h>
+#include <numerics/mesh_worker_info.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace LocalIntegrators
+{
+/**
+ * @brief Local integrators related to elasticity problems.
+ *
+ * @ingroup Integrators
+ * @author Guido Kanschat
+ * @date 2010
+ */
+  namespace Elasticity
+  {
+                                    /**
+                                     * The weak form of the grad-div
+                                     * operator penalizing volume changes
+                                     * @f[
+                                     * \int_Z \nabla\!\cdot\!u
+                                     * \nabla\!\cdot\!v \,dx
+                                     * @f]
+                                     */
+    template <int dim>
+    void grad_div_matrix (
+      FullMatrix<double>& M,
+      const FEValuesBase<dim>& fe,
+      double factor = 1.)
+    {
+      const unsigned int n_dofs = fe.dofs_per_cell;
+      
+      AssertDimension(fe.get_fe().n_components(), dim);
+      AssertDimension(M.m(), n_dofs);
+      AssertDimension(M.n(), n_dofs);
+      
+      for (unsigned k=0;k<fe.n_quadrature_points;++k)
+       {
+         const double dx = factor * fe.JxW(k);
+         for (unsigned i=0;i<n_dofs;++i)
+           for (unsigned j=0;j<n_dofs;++j)
+             {
+               double dv = 0.;
+               double du = 0.;
+               for (unsigned int d=0;d<dim;++d)
+                 {
+                   dv += fe.shape_grad_component(i,k,d)[d];
+                   du += fe.shape_grad_component(j,k,d)[d];
+                 }
+               
+               M(i,j) += dx * du * dv;
+             }
+       }
+    }
+  }
+}
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 37f32d258c41ae3c7a759022168e732755e2fca8..6b3bc2fa3f63db5c1b48e5fc852817d9ad3757a4 100644 (file)
@@ -29,6 +29,7 @@ namespace LocalIntegrators
 /**
  * @brief Local integrators related to <i>L<sup>2</sup></i>-inner products.
  *
+ * @ingroup Integrators
  * @author Guido Kanschat
  * @date 2010
  */
index 88a873f1f2b0f22bdd439b0fe4c3803f5a899ee1..9d213be76d4204ccf5565476251806f2ffbcbc32 100644 (file)
@@ -44,6 +44,7 @@ namespace LocalIntegrators
  * the latter case, the Laplacian is applied to each component
  * separately.
  *
+ * @ingroup Integrators
  * @author Guido Kanschat
  * @date 2008, 2009, 2010
  */
index bbd5c1881b37760186787a0feffe5f0d6da38b92..f7ed4a2eb7f7dff46b2aa36454760aabc6690eb8 100644 (file)
@@ -118,6 +118,8 @@ void MatrixIntegrator<dim>::cell(
 }
  * @endcode
  * See step-42 for a worked out example of this code.
+ *
+ * @ingroup Integrators
  */
 namespace LocalIntegrators
 {
diff --git a/deal.II/include/deal.II/integrators/maxwell.h b/deal.II/include/deal.II/integrators/maxwell.h
new file mode 100644 (file)
index 0000000..8b665fa
--- /dev/null
@@ -0,0 +1,135 @@
+//---------------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//---------------------------------------------------------------------------
+#ifndef __deal2__integrators_maxwell_h
+#define __deal2__integrators_maxwell_h
+
+
+#include <base/config.h>
+#include <base/exceptions.h>
+#include <base/quadrature.h>
+#include <lac/full_matrix.h>
+#include <fe/mapping.h>
+#include <fe/fe_values.h>
+#include <numerics/mesh_worker_info.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace LocalIntegrators
+{
+/**
+ * @brief Local integrators related to curl operators and their traces.
+ *
+ * @ingroup Integrators
+ * @author Guido Kanschat
+ * @date 2010
+ */
+  namespace Maxwell
+  {
+                                    /**
+                                     * The curl-curl operator
+                                     * @f[
+                                     * \int_Z \nabla\!\times\! u \cdot
+                                     * \nabla\!\times\! v \,dx
+                                     * @f]
+                                     */
+    template <int dim>
+    void curl_curl_matrix (
+      FullMatrix<double>& M,
+      const FEValuesBase<dim>& fe,
+      const double factor = 1.)
+    {
+      const unsigned int n_dofs = fe.dofs_per_cell;
+      
+      AssertDimension(fe.get_fe().n_components(), dim);
+      AssertDimension(M.m(), n_dofs);
+      AssertDimension(M.n(), n_dofs);
+      
+                                      // Depending on the dimension,
+                                      // the cross product is either
+                                      // a scalar (2d) or a vector
+                                      // (3d). Accordingly, in the
+                                      // latter case we have to sum
+                                      // up three bilinear forms, but
+                                      // in 2d, we don't. Thus, we
+                                      // need to adapt the loop over
+                                      // all dimensions
+      const unsigned int d_max = (dim==2) ? 1 : dim;
+      
+      for (unsigned k=0;k<fe.n_quadrature_points;++k)
+       {
+         const double dx = factor * fe.JxW(k);
+         for (unsigned i=0;i<n_dofs;++i)
+           for (unsigned j=0;j<n_dofs;++j)
+             for (unsigned int d=0;d<d_max;++d)
+               {
+                 const unsigned int d1 = (d+1)%dim;
+                 const unsigned int d2 = (d+2)%dim;
+                 
+                 const double cv = fe.shape_grad_component(i,k,d1)[d2] - fe.shape_grad_component(i,k,d2)[d1];
+                 const double cu = fe.shape_grad_component(j,k,d1)[d2] - fe.shape_grad_component(j,k,d2)[d1];
+                 
+                 M(i,j) += dx * cu * cv;
+               }
+       }
+    }
+    
+                                    /**
+                                     * The curl operator
+                                     * @f[
+                                     * \int_Z \nabla\!\times\! u \cdot v \,dx.
+                                     * @f]
+                                     *
+                                     * This is the standard curl
+                                     * operator in 3D and the scalar
+                                     * curl in 2D.
+                                     */
+    template <int dim>
+    void curl_matrix (
+      FullMatrix<double>& M,
+      const FEValuesBase<dim>& fe,
+      const FEValuesBase<dim>& fetest,
+      double factor = 1.)
+    {
+      unsigned int t_comp = (dim==3) ? dim : 1;
+      const unsigned int n_dofs = fe.dofs_per_cell;
+      const unsigned int t_dofs = fetest.dofs_per_cell;
+      AssertDimension(fe.get_fe().n_components(), dim);
+      AssertDimension(fetest.get_fe().n_components(), t_comp);
+      AssertDimension(M.m(), t_dofs);
+      AssertDimension(M.n(), n_dofs);
+      
+      const unsigned int d_max = (dim==2) ? 1 : dim;
+      
+      for (unsigned k=0;k<fe.n_quadrature_points;++k)
+       {
+         const double dx = fe.JxW(k) * factor;
+         for (unsigned i=0;i<t_dofs;++i)
+           for (unsigned j=0;j<n_dofs;++j)
+             for (unsigned int d=0;d<d_max;++d)
+               {
+                 const unsigned int d1 = (d+1)%dim;
+                 const unsigned int d2 = (d+2)%dim;
+                 
+                 const double vv = fetest.shape_value_component(i,k,d);
+                 const double cu = fe.shape_grad_component(j,k,d1)[d2] - fe.shape_grad_component(j,k,d2)[d1];
+                 M(i,j) += dx * cu * vv;
+               }
+       }
+    }
+  }
+}
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 2321eb8f58efb171b43601d79f56c6d74a4c6e9b..8ea5441e97e02abbd5cd60af76c5fc11f362ab69 100644 (file)
@@ -167,7 +167,9 @@ template<int,int> class MGDoFHandler;
  * we may need more than the default update flags.
  *
  * @ingroup MeshWorker
- * @author Guido Kanschat, 2009
+ * @ingroup Integrators
+ * @author Guido Kanschat
+ * @date 2009
  */
 namespace MeshWorker
 {

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.