00001 /** \page FEDoc Finite elements concepts in Sundance 00002 00003 <BODY> 00004 <H1> Finite elements concepts in Sundance </H1> 00005 00006 00007 In the finite-element method, a PDE is 00008 written in weak form on an infinite-dimensional space 00009 and then discretized onto 00010 a finite-dimensional space as specified by the choice of Mesh and 00011 by the choice of basis function. 00012 00013 There are many possible weak formulations for any given PDE, for instance 00014 Galerkin or first-order least squares. Furthermore, there is then 00015 some freedom in choice of basis function (often constrained by 00016 stability considerations). Sundance is designed to allow great flexibility 00017 in choice of formulation and discretization method. 00018 00019 \section FunctionSpace Function spaces and basis families 00020 00021 The unknown and variational (test) functions in a weak equation are 00022 represented by UnknownFunction and TestFunction objects, 00023 respectively. The first constructor argument to an UnknownFunction 00024 or TestFunction object is the BasisFamily with which it 00025 will be discretized. For example, 00026 \code 00027 Expr T = new UnknownFunction(new Lagrange(1), "T"); 00028 \endcode 00029 creates a function T that will be discretized with first-order Lagrange 00030 interpolation. 00031 00032 A BasisFamily such as Lagrange 00033 is a geometry-independent specification of a basis function; the 00034 basis family object 00035 automatically hook to the low-level functions appropriate to 00036 the cell being integrated (e.g. \f$ x \f$ and \f$1-x\f$ for a line, and 00037 so on). 00038 00039 00040 00041 \section VariationalForms Weak equations 00042 00043 The continuous formulation of a finite-elements problem is as a 00044 weak equation, for example 00045 \f[ - \int_\Omega \nabla u \cdot \nabla v + \int_{\Gamma_1} v G 00046 - \int_\Omega v f = 0 \f] 00047 plus possible restrictions on the variations on specified 00048 boundary segments (see \ref BoundaryConditions below). 00049 00050 To represent a weak equation, use the Sundance Integral object. Assuming 00051 that the boundary segment \f$\Gamma_1\f$ has been represented 00052 with a CellSet called gamma 00053 the equation above would be written 00054 \code 00055 Integral myEqn = Integral(-(grad*v)*(grad*u)) 00056 + Integral(-v*f) 00057 + Integral(gamma, v*G); 00058 \endcode 00059 Integrals with no domain specified are taken to be over all maximal cells. 00060 Notice that Integral objects can be added using the overloaded \c+\c operator. 00061 00062 Since the integrals will be done with quadrature, you can give a specification 00063 of a quadrature family as an optional argument to Integral. For example, 00064 \code 00065 Integral myEqn = Integral(someCellSet, v*sin(x), new GaussLegendre(8)); 00066 \endcode 00067 will cause 8th-order Gauss-Legendre quadrature to be used on whatever 00068 cells get integrated. You can use this feature on a term-by-term 00069 basis to assign high-order 00070 quadrature for complicated terms but low-order quadrature for 00071 all others, for example 00072 \code 00073 Integral myEqn = Integral(v*sin(4*x), new GaussLegendre(16)) 00074 + Integral(v*u, new GaussLegendre(2)); 00075 \endcode 00076 The default quadrature rule is 4th-order Gauss-Legendre. 00077 00078 \warning 00079 The Sundance Integral object is essentially just a bookkeeper that groups 00080 domains with expressions. It 00081 doesn't know how to do integration by parts, so 00082 for problems in which integration by parts is necessary 00083 you must do that step by hand before putting your problem into an Integral 00084 object. 00085 00086 00087 \section BoundaryConditions Boundary conditions 00088 00089 There are many ways to handle boundary conditions in finite-elements 00090 problems. 00091 00092 In Galerkin methods, Neumann BCs are easy: the BC 00093 \f${\hat n} \cdot \nabla u = G \f$ along a CellSet top 00094 is incorporated by adding 00095 \code 00096 Integral(top, v*G) 00097 \endcode 00098 to the weak equation, where v is the variation of u and G is the 00099 expression which represents the flux. Natural BCs, which are zero-flux 00100 Neumann BCs, are especially easy, since the integral is indentically 00101 zero. To do natural BCs, you do nothing at all. 00102 00103 Dirichlet BCs can done by replacing the weak equation on the boundary 00104 nodes with another weak equation representing the BCs. In Sundance, this 00105 is done with the EssentialBC object. To apply the weak condition 00106 \f$u=sin(y)\f$ on the CellSet left, write 00107 \code 00108 EssentialBC bc = EssentialBC(left, v*(u-sin(y))); 00109 \endcode 00110 EssentialBC objects can be combined using the "and" operator 00111 \&\&, for example 00112 \code 00113 EssentialBC bc = EssentialBC(left, v*(u-sin(y))) 00114 && EssentialBC(right, v*(u-cos(y))); 00115 \endcode 00116 00117 In least-squares formulations, you can do dirichlet BCs in least-squares 00118 form (i.e., as just another Integral) or in essential form. Using 00119 the least-squares form has the advantage that your matrix remains SPD, but 00120 the BCs will not always be satisfied exactly. 00121 00122 You can also 00123 do Dirichlet BCs by introducing lagrange multipliers along the boundary. 00124 In this case, the BC becomes just another integral to be added to 00125 the weak equation. 00126 This method is sometimes useful, but produces indefinite matrices. 00127 00128 00129 <H2> 00130 Next section: \ref LinearSolvers 00131 </h2> 00132 00133 </BODY> 00134 */ 00135 00136 00137 00138 00139 00140