From: kanschat Date: Tue, 9 Nov 2010 15:00:06 +0000 (+0000) Subject: add local integrators and list namespaces in module X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=dce20a6af71b2f329671e521312fbcf3d0ad8d73;p=dealii-svn.git add local integrators and list namespaces in module git-svn-id: https://svn.dealii.org/trunk@22642 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/include/deal.II/integrators/differential.h b/deal.II/include/deal.II/integrators/differential.h index 5beec751cb..00af25a290 100644 --- a/deal.II/include/deal.II/integrators/differential.h +++ b/deal.II/include/deal.II/integrators/differential.h @@ -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 index 0000000000..f30d96558d --- /dev/null +++ b/deal.II/include/deal.II/integrators/elasticity.h @@ -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 +#include +#include +#include +#include +#include +#include + +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 + void grad_div_matrix ( + FullMatrix& M, + const FEValuesBase& 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;kL2-inner products. * + * @ingroup Integrators * @author Guido Kanschat * @date 2010 */ diff --git a/deal.II/include/deal.II/integrators/laplace.h b/deal.II/include/deal.II/integrators/laplace.h index 88a873f1f2..9d213be76d 100644 --- a/deal.II/include/deal.II/integrators/laplace.h +++ b/deal.II/include/deal.II/integrators/laplace.h @@ -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 */ diff --git a/deal.II/include/deal.II/integrators/local_integrators.h b/deal.II/include/deal.II/integrators/local_integrators.h index bbd5c1881b..f7ed4a2eb7 100644 --- a/deal.II/include/deal.II/integrators/local_integrators.h +++ b/deal.II/include/deal.II/integrators/local_integrators.h @@ -118,6 +118,8 @@ void MatrixIntegrator::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 index 0000000000..8b665fa2fc --- /dev/null +++ b/deal.II/include/deal.II/integrators/maxwell.h @@ -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 +#include +#include +#include +#include +#include +#include + +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 + void curl_curl_matrix ( + FullMatrix& M, + const FEValuesBase& 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 + void curl_matrix ( + FullMatrix& M, + const FEValuesBase& fe, + const FEValuesBase& 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 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 {