Line data Source code
1 : /*
2 : * Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Christian Wehner
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__CROUZEIX_RAVIART__
34 : #define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__CROUZEIX_RAVIART__
35 :
36 : #include "lib_grid/grid/grid_base_objects.h"
37 : #include "lib_disc/reference_element/reference_element_util.h"
38 : #include "lib_disc/local_finite_element/local_shape_function_set.h"
39 :
40 : namespace ug{
41 :
42 :
43 : /// Crouzeix - Raviart Set
44 : template <typename TRefElem>
45 : class CrouzeixRaviartBase
46 : {
47 : public:
48 : /// dimension of reference element
49 : static const int dim = TRefElem::dim;
50 :
51 : /// Number of shape functions
52 : static const size_t nsh = TRefElem::numSides;
53 :
54 : public:
55 : /// constructor
56 0 : CrouzeixRaviartBase()
57 0 : {
58 0 : for(size_t i = 0; i < nsh; ++i)
59 0 : m_vLocalDoF[i] = LocalDoF(dim-1, i, 0);
60 : }
61 :
62 : /// returns the type of reference element
63 : ReferenceObjectID roid() const {return TRefElem::REFERENCE_OBJECT_ID;}
64 :
65 : /// returns the total number of DoFs on the finite element
66 : size_t num_dof() const {return nsh;};
67 : size_t num_sh() const {return nsh;}
68 :
69 : /// returns the number of DoFs on a sub-geometric object type
70 : size_t num_dof(ReferenceObjectID type) const
71 : {
72 0 : if(ReferenceElementDimension(type) == dim-1) return 1;
73 0 : else return 0;
74 : }
75 :
76 : /// returns the dof storage
77 0 : const LocalDoF& local_dof(size_t dof) const {return m_vLocalDoF[dof];}
78 :
79 : /// returns if the local dof position are exact
80 : bool exact_position_available() const {return true;};
81 :
82 : /// \copydoc ug::LocalShapeFunctionSet::continuous()
83 : bool continuous() const {return false;}
84 :
85 : protected:
86 : /// association to elements
87 : LocalDoF m_vLocalDoF[nsh];
88 : };
89 :
90 : /// Lagrange Shape Function Set without virtual functions and fixed order
91 : template <typename TRefElem>
92 : class CrouzeixRaviartLSFS;
93 :
94 : ////////////////////////////////////////////////////////////////////////////////
95 : // ReferenceTriangle
96 : //
97 : // function space span {1,x,y}
98 : //
99 : ////////////////////////////////////////////////////////////////////////////////
100 :
101 : template <>
102 : class CrouzeixRaviartLSFS<ReferenceTriangle>
103 : : public CrouzeixRaviartBase<ReferenceTriangle>,
104 : public BaseLSFS<CrouzeixRaviartLSFS<ReferenceTriangle>, 2>
105 : {
106 : public:
107 : /// Dimension, where shape functions are defined
108 : static const int dim = 2;
109 :
110 : public:
111 : /// \copydoc ug::LocalShapeFunctionSet::position()
112 0 : bool position(size_t i, MathVector<dim>& pos) const
113 : {
114 0 : switch(i)
115 : {
116 0 : case 0: pos[0] = 0.5;
117 0 : pos[1] = 0.0; return true;
118 0 : case 1: pos[0] = 0.5;
119 0 : pos[1] = 0.5; return true;
120 0 : case 2: pos[0] = 0.0;
121 0 : pos[1] = 0.5; return true;
122 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
123 : " not found. Only "<<nsh<<" shapes present.");
124 : }
125 : }
126 :
127 : /// \copydoc ug::LocalShapeFunctionSet::shape()
128 0 : number shape(const size_t i, const MathVector<dim>& x) const
129 : {
130 0 : switch(i)
131 : {
132 0 : case 0: return 1-2*x[1];
133 0 : case 1: return -1+2*x[0]+2*x[1];
134 0 : case 2: return 1-2*x[0];
135 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
136 : " not found. Only "<<nsh<<" shapes present.");
137 : }
138 : }
139 :
140 : /// \copydoc ug::LocalShapeFunctionSet::grad()
141 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
142 : {
143 0 : switch(i)
144 : {
145 0 : case 0: g[0] = 0.0;
146 0 : g[1] = -2.0; return;
147 0 : case 1: g[0] = 2.0;
148 0 : g[1] = 2.0; return;
149 0 : case 2: g[0] = -2.0;
150 0 : g[1] = 0.0; return;
151 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
152 : " not found. Only "<<nsh<<" shapes present.");
153 : }
154 : }
155 : };
156 :
157 :
158 : ////////////////////////////////////////////////////////////////////////////////
159 : // ReferenceQuadrilateral
160 : //
161 : // function space span {1,x,y,x^2-y^2}
162 : //
163 : ////////////////////////////////////////////////////////////////////////////////
164 :
165 : template <>
166 : class CrouzeixRaviartLSFS<ReferenceQuadrilateral>
167 : : public CrouzeixRaviartBase<ReferenceQuadrilateral>,
168 : public BaseLSFS<CrouzeixRaviartLSFS<ReferenceQuadrilateral>, 2>
169 : {
170 : public:
171 : /// Dimension, where shape functions are defined
172 : static const int dim = 2;
173 :
174 : public:
175 : /// \copydoc ug::LocalShapeFunctionSet::position()
176 0 : bool position(size_t i, MathVector<dim>& pos) const
177 : {
178 0 : switch(i)
179 : {
180 0 : case 0: pos[0] = 0.5;
181 0 : pos[1] = 0.0; return true;
182 0 : case 1: pos[0] = 1.0;
183 0 : pos[1] = 0.5; return true;
184 0 : case 2: pos[0] = 0.5;
185 0 : pos[1] = 1.0; return true;
186 0 : case 3: pos[0] = 0.0;
187 0 : pos[1] = 0.5; return true;
188 : return true;
189 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
190 : " not found. Only "<<nsh<<" shapes present.");
191 : }
192 : }
193 :
194 : /// \copydoc ug::LocalShapeFunctionSet::shape()
195 0 : number shape(const size_t i, const MathVector<dim>& x) const
196 : {
197 0 : switch(i)
198 : {
199 0 : case 0: return 0.75+x[0]-2*x[1]-x[0]*x[0]+x[1]*x[1];
200 0 : case 1: return -0.25+x[1]+x[0]*x[0]-x[1]*x[1];
201 0 : case 2: return -0.25+x[0]-x[0]*x[0]+x[1]*x[1];
202 0 : case 3: return 0.75-2*x[0]+x[1]+x[0]*x[0]-x[1]*x[1];
203 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
204 : " not found. Only "<<nsh<<" shapes present.");
205 : }
206 : }
207 :
208 : /// \copydoc ug::LocalShapeFunctionSet::grad()
209 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
210 : {
211 0 : switch(i)
212 : {
213 0 : case 0: g[0] = 1-2*x[0];
214 0 : g[1] = -2.0+2*x[1]; return;
215 0 : case 1: g[0] = 2.0*x[0];
216 0 : g[1] = 1-2.0*x[1]; return;
217 0 : case 2: g[0] = 1-2.0*x[0];
218 0 : g[1] = 2*x[1]; return;
219 0 : case 3: g[0] = -2+2*x[0];
220 0 : g[1] = 1-2*x[1]; return;
221 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
222 : " not found. Only "<<nsh<<" shapes present.");
223 : }
224 : }
225 : };
226 :
227 : ////////////////////////////////////////////////////////////////////////////////
228 : // ReferenceTetrahedron
229 : //
230 : // function space span {1,x,y,z}
231 : //
232 : ////////////////////////////////////////////////////////////////////////////////
233 :
234 : template <>
235 : class CrouzeixRaviartLSFS<ReferenceTetrahedron>
236 : : public CrouzeixRaviartBase<ReferenceTetrahedron>,
237 : public BaseLSFS<CrouzeixRaviartLSFS<ReferenceTetrahedron>, 3>
238 : {
239 : public:
240 : /// Dimension, where shape functions are defined
241 : static const int dim = ReferenceTetrahedron::dim;
242 :
243 : public:
244 : /// \copydoc ug::LocalShapeFunctionSet::position()
245 0 : bool position(size_t i, MathVector<dim>& pos) const
246 : {
247 0 : switch(i)
248 : {
249 0 : case 0: pos[0] = 1.0/3.0;
250 0 : pos[1] = 1.0/3.0;
251 0 : pos[2] = 0.0; return true;
252 0 : case 1: pos[0] = 1.0/3.0;
253 0 : pos[1] = 1.0/3.0;
254 0 : pos[2] = 1.0/3.0; return true;
255 0 : case 2: pos[0] = 0.0;
256 0 : pos[1] = 1.0/3.0;
257 0 : pos[2] = 1.0/3.0; return true;
258 0 : case 3: pos[0] = 1.0/3.0;
259 0 : pos[1] = 0.0;
260 0 : pos[2] = 1.0/3.0; return true;
261 : return true;
262 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
263 : " not found. Only "<<nsh<<" shapes present.");
264 : }
265 : }
266 :
267 : /// \copydoc ug::LocalShapeFunctionSet::shape()
268 0 : number shape(const size_t i, const MathVector<dim>& x) const
269 : {
270 0 : switch(i)
271 : {
272 0 : case 0: return 1 - 3*x[2];;
273 0 : case 1: return -2 + 3*x[0]+3*x[1]+3*x[2];
274 0 : case 2: return 1 - 3*x[0];
275 0 : case 3: return 1 - 3*x[1];
276 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
277 : " not found. Only "<<nsh<<" shapes present.");
278 : }
279 : }
280 :
281 : /// \copydoc ug::LocalShapeFunctionSet::grad()
282 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
283 : {
284 0 : switch(i)
285 : {
286 0 : case 0: g[0] = 0.0;
287 0 : g[1] = 0.0;
288 0 : g[2] = -3.0;return;
289 0 : case 1: g[0] = 3.0;
290 0 : g[1] = 3.0;
291 0 : g[2] = 3.0;return;
292 0 : case 2: g[0] = -3.0;
293 0 : g[1] = 0.0;
294 0 : g[2] = 0.0;return;
295 0 : case 3: g[0] = 0.0;
296 0 : g[1] = -3;
297 0 : g[2] = 0.0;return;
298 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
299 : " not found. Only "<<nsh<<" shapes present.");
300 : }
301 : }
302 : };
303 :
304 :
305 : ////////////////////////////////////////////////////////////////////////////////
306 : // ReferenceHexahedron
307 : //
308 : // function space span {1,x,y,z,x^2-y^2,y^2-z^2}
309 : //
310 : ////////////////////////////////////////////////////////////////////////////////
311 :
312 : template <>
313 : class CrouzeixRaviartLSFS<ReferenceHexahedron>
314 : : public CrouzeixRaviartBase<ReferenceHexahedron>,
315 : public BaseLSFS<CrouzeixRaviartLSFS<ReferenceHexahedron>, 3>
316 : {
317 : public:
318 : /// Dimension, where shape functions are defined
319 : static const int dim = ReferenceHexahedron::dim;
320 :
321 : public:
322 : /// \copydoc ug::LocalShapeFunctionSet::position()
323 0 : bool position(size_t i, MathVector<dim>& pos) const
324 : {
325 0 : switch(i)
326 : {
327 0 : case 0: pos[0] = 0.5;
328 0 : pos[1] = 0.5;
329 0 : pos[2] = 0.0;return true;
330 0 : case 1: pos[0] = 0.5;
331 0 : pos[1] = 0.0;
332 0 : pos[2] = 0.5; return true;
333 0 : case 2: pos[0] = 1.0;
334 0 : pos[1] = 0.5;
335 0 : pos[2] = 0.5; return true;
336 0 : case 3: pos[0] = 0.5;
337 0 : pos[1] = 1.0;
338 0 : pos[2] = 0.5; return true;
339 0 : case 4: pos[0] = 0.0;
340 0 : pos[1] = 0.5;
341 0 : pos[2] = 0.5; return true;
342 0 : case 5: pos[0] = 0.5;
343 0 : pos[1] = 0.5;
344 0 : pos[2] = 1.0; return true;
345 : return true;
346 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
347 : " not found. Only "<<nsh<<" shapes present.");
348 : }
349 : }
350 :
351 : /// \copydoc ug::LocalShapeFunctionSet::shape()
352 0 : number shape(const size_t i, const MathVector<dim>& x) const
353 : {
354 0 : switch(i)
355 : {
356 0 : case 0: return 1.0/3.0*(2 + 2*x[0] + 2*x[1] - 7*x[2]-2*x[0]*x[0]-2*x[1]*x[1]+4*x[2]*x[2]);
357 0 : case 1: return 1.0/3.0*(2+2*x[0]-7*x[1]+2*x[2]-2*x[0]*x[0]+4*x[1]*x[1]-2*x[2]*x[2]);
358 0 : case 2: return 1.0/3.0*(-1-x[0]+2*x[1]+2*x[2]+4*x[0]*x[0]-2*x[1]*x[1]-2*x[2]*x[2]);
359 0 : case 3: return 1.0/3.0*(-1+2*x[0]-x[1]+2*x[2]-2*x[0]*x[0]+4*x[1]*x[1]-2*x[2]*x[2]);
360 0 : case 4: return 1.0/3.0*(2-7*x[0]+2*x[1]+2*x[2]+4*x[0]*x[0]-2*x[1]*x[1]-2*x[2]*x[2]);
361 0 : case 5: return 1.0/3.0*(-1+2*x[0]+2*x[1]-x[2]-2*x[0]*x[0]-2*x[1]*x[1]+4*x[2]*x[2]);
362 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
363 : " not found. Only "<<nsh<<" shapes present.");
364 : }
365 : }
366 :
367 : /// \copydoc ug::LocalShapeFunctionSet::grad()
368 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
369 : {
370 0 : switch(i)
371 : {
372 0 : case 0: g[0] = 1.0/3.0*(2-4*x[0]);
373 0 : g[1] = 1.0/3.0*(2-4*x[1]);
374 0 : g[2] = 1.0/3.0*(-7+8*x[2]);return;
375 0 : case 1: g[0] = 1.0/3.0*(2-4*x[0]);
376 0 : g[1] = 1.0/3.0*(-7+8*x[1]);
377 0 : g[2] = 1.0/3.0*(2-4*x[2]);return;
378 0 : case 2: g[0] = 1.0/3.0*(-1+8*x[0]);
379 0 : g[1] = 1.0/3.0*(2-4*x[1]);
380 0 : g[2] = 1.0/3.0*(2-4*x[2]);return;
381 0 : case 3: g[0] = 1.0/3.0*(2-4*x[0]);
382 0 : g[1] = 1.0/3.0*(-1+8*x[1]);
383 0 : g[2] = 1.0/3.0*(2-4*x[2]);return;
384 0 : case 4: g[0] = 1.0/3.0*(-7+8*x[0]);
385 0 : g[1] = 1.0/3.0*(2-4*x[1]);
386 0 : g[2] = 1.0/3.0*(2-4*x[2]);return;
387 0 : case 5: g[0] = 1.0/3.0*(2-4*x[0]);
388 0 : g[1] = 1.0/3.0*(2-4*x[1]);
389 0 : g[2] = 1.0/3.0*(-1+8*x[2]);return;
390 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
391 : " not found. Only "<<nsh<<" shapes present.");
392 : }
393 : }
394 : };
395 :
396 : ////////////////////////////////////////////////////////////////////////////////
397 : // ReferencePrism
398 : //
399 : // function space span {1,x,y,z,x^2-y^2-z^2}
400 : //
401 : ////////////////////////////////////////////////////////////////////////////////
402 :
403 : template <>
404 : class CrouzeixRaviartLSFS<ReferencePrism>
405 : : public CrouzeixRaviartBase<ReferencePrism>,
406 : public BaseLSFS<CrouzeixRaviartLSFS<ReferencePrism>, 3>
407 : {
408 : public:
409 : /// Dimension, where shape functions are defined
410 : static const int dim = 3;
411 :
412 : public:
413 : /// \copydoc ug::LocalShapeFunctionSet::position()
414 0 : bool position(size_t i, MathVector<dim>& pos) const
415 : {
416 0 : switch(i)
417 : {
418 0 : case 0: pos[0] = 1.0/3.0;
419 0 : pos[1] = 1.0/3.0;
420 0 : pos[2] = 0.0;return true;
421 0 : case 1: pos[0] = 0.5;
422 0 : pos[1] = 0.0;
423 0 : pos[2] = 0.5; return true;
424 0 : case 2: pos[0] = 0.5;
425 0 : pos[1] = 0.5;
426 0 : pos[2] = 0.5; return true;
427 0 : case 3: pos[0] = 0.0;
428 0 : pos[1] = 0.5;
429 0 : pos[2] = 0.5; return true;
430 0 : case 4: pos[0] = 1.0/3.0;
431 0 : pos[1] = 1.0/3.0;
432 0 : pos[2] = 1.0; return true;
433 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
434 : " not found. Only "<<nsh<<" shapes present.");
435 : }
436 : }
437 :
438 : /// \copydoc ug::LocalShapeFunctionSet::shape()
439 0 : number shape(const size_t i, const MathVector<dim>& x) const
440 : {
441 0 : switch(i)
442 : {
443 0 : case 0: return 1+x[0]-x[1]-3*x[2]-2*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
444 0 : case 1: return 1.0/3.0*(2-2*x[0]-4*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
445 0 : case 2: return 1.0/3.0*(-4+4*x[0]+8*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
446 0 : case 3: return 1.0/3.0*(2-8*x[0]+2*x[1]+4*x[2]+4*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]));
447 0 : case 4: return x[0]-x[1]-x[2]-2*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
448 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
449 : " not found. Only "<<nsh<<" shapes present.");
450 : }
451 : }
452 :
453 : /// \copydoc ug::LocalShapeFunctionSet::grad()
454 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
455 : {
456 0 : switch(i)
457 : {
458 0 : case 0: g[0] = 1-4*x[0];
459 0 : g[1] = -1+4*x[1];
460 0 : g[2] = -3+4*x[2];return;
461 0 : case 1: g[0] = 1.0/3.0*(-2+8*x[0]);
462 0 : g[1] = 1.0/3.0*(-4-8*x[1]);
463 0 : g[2] = 1.0/3.0*(4-8*x[2]);return;
464 0 : case 2: g[0] = 1.0/3.0*(4+8*x[0]);
465 0 : g[1] = 1.0/3.0*(8-8*x[1]);
466 0 : g[2] = 1.0/3.0*(4-8*x[2]);return;
467 0 : case 3: g[0] = 1.0/3.0*(-8+8*x[0]);
468 0 : g[1] = 1.0/3.0*(2-8*x[1]);
469 0 : g[2] = 1.0/3.0*(4-8*x[2]);return;
470 0 : case 4: g[0] = 1-4*x[0];
471 0 : g[1] = -1+4*x[1];
472 0 : g[2] = -1+4*x[2];return;
473 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
474 : " not found. Only "<<nsh<<" shapes present.");
475 : }
476 : }
477 : };
478 :
479 :
480 : ////////////////////////////////////////////////////////////////////////////////
481 : // ReferencePyramid
482 : //
483 : // function space span {1,x,y,z,x^2-y^2-z^2}
484 : //
485 : ////////////////////////////////////////////////////////////////////////////////
486 :
487 : template <>
488 : class CrouzeixRaviartLSFS<ReferencePyramid>
489 : : public CrouzeixRaviartBase<ReferencePyramid>,
490 : public BaseLSFS<CrouzeixRaviartLSFS<ReferencePyramid>, 3>
491 : {
492 : public:
493 : /// Dimension, where shape functions are defined
494 : static const int dim = 3;
495 :
496 : public:
497 : /// \copydoc ug::LocalShapeFunctionSet::position()
498 0 : bool position(size_t i, MathVector<dim>& pos) const
499 : {
500 0 : switch(i)
501 : {
502 0 : case 0: pos[0] = 0.5;
503 0 : pos[1] = 0.5;
504 0 : pos[2] = 0.0;return true;
505 0 : case 1: pos[0] = 1.0/3.0;
506 0 : pos[1] = 0.0;
507 0 : pos[2] = 1.0/3.0; return true;
508 0 : case 2: pos[0] = 2.0/3.0;
509 0 : pos[1] = 1.0/3.0;
510 0 : pos[2] = 1.0/3.0; return true;
511 0 : case 3: pos[0] = 1.0/3.0;
512 0 : pos[1] = 2.0/3.0;
513 0 : pos[2] = 1.0/3.0; return true;
514 0 : case 4: pos[0] = 0.0;
515 0 : pos[1] = 1.0/3.0;
516 0 : pos[2] = 1.0/3.0; return true;
517 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
518 : " not found. Only "<<nsh<<" shapes present.");
519 : }
520 : }
521 :
522 : /// \copydoc ug::LocalShapeFunctionSet::shape()
523 0 : number shape(const size_t i, const MathVector<dim>& x) const
524 : {
525 0 : switch(i)
526 : {
527 0 : case 0: return 1-3*x[2];
528 0 : case 1: return 0.75+1.5*x[0]-3*x[1]-0.75*x[2]-2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
529 0 : case 2: return -0.75+1.5*x[1]+2.25*x[2]+2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
530 0 : case 3: return -0.75+1.5*x[0]+0.75*x[2]-2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
531 0 : case 4: return 0.75-3*x[0]+1.5*x[1]+0.75*x[2]+2.25*(x[0]*x[0]-x[1]*x[1]-x[2]*x[2]);
532 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
533 : " not found. Only "<<nsh<<" shapes present.");
534 : }
535 : }
536 :
537 : /// \copydoc ug::LocalShapeFunctionSet::grad()
538 0 : void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
539 : {
540 0 : switch(i)
541 : {
542 0 : case 0: g[0] = 0.0;
543 0 : g[1] = 0.0;
544 0 : g[2] = -3.0;return;
545 0 : case 1: g[0] = 1.5-4.5*x[0];
546 0 : g[1] = -3+4.5*x[1];
547 0 : g[2] = -0.75+4.5*x[2];return;
548 0 : case 2: g[0] = 4.5*x[0];
549 0 : g[1] = 1.5-4.5*x[1];
550 0 : g[2] = 2.25-4.5*x[2];return;
551 0 : case 3: g[0] = 1.5-4.5*x[0];
552 0 : g[1] = 4.5*x[1];
553 0 : g[2] = 0.75+4.5*x[2];return;
554 0 : case 4: g[0] = -3+4.5*x[0];
555 0 : g[1] = 1.5-4.5*x[1];
556 0 : g[2] = 0.75-4.5*x[2];return;
557 0 : default: UG_THROW("CrouzeixRaviartLSFS: shape function "<<i<<
558 : " not found. Only "<<nsh<<" shapes present.");
559 : }
560 : }
561 : };
562 :
563 :
564 : } //namespace ug
565 :
566 : #endif /* __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__CROUZEIX_RAVIART__ */
567 :
|