src/fem/nedelecelement.h

Go to the documentation of this file.
00001 /***************************************************************************
00002                           nedelecelement.h  -  description
00003                              -------------------
00004     begin                : Fri Dec 19 2003
00005     copyright            : (C) 2003 by Roman Geus
00006     email                : roman.geus@psi.ch
00007 ***************************************************************************/
00008 
00009 /***************************************************************************
00010  *                                                                         *
00011  *   This program is free software; you can redistribute it and/or modify  *
00012  *   it under the terms of the GNU General Public License as published by  *
00013  *   the Free Software Foundation; either version 2 of the License, or     *
00014  *   (at your option) any later version.                                   *
00015  *                                                                         *
00016  ***************************************************************************/
00017 
00018 #ifndef NEDELECELEMENT_H
00019 #define NEDELECELEMENT_H
00020 
00021 #include "tetmesh/tet.h"
00022 #include "tetmesh/tetmesh.h"
00023 #include "tetmesh/vector4.h"
00024 #include "tetmesh/vector3.h"
00025 #include "colarray/vector.h"
00026 #include "colarray/matrix.h"
00027 
00028 using namespace colarray;
00029 
00037 class NedelecElement {
00038 public: 
00039     NedelecElement();
00040     virtual ~NedelecElement();
00042     virtual void get_element_matrix(mesh::Tet* tet, colarray::Matrix<double>& Ae, colarray::Matrix<double>& Me) = 0;
00044     virtual mesh::Vector3 eval(mesh::Tet* tet, const mesh::Vector3& x, const colarray::Vector<double>& dof) = 0;
00046     virtual mesh::Vector3 eval_curl(mesh::Tet* tet, const mesh::Vector3& x, const colarray::Vector<double>& dof) = 0;
00048     virtual int get_nof_dofs() const = 0;
00050     virtual int get_order() const = 0;
00054     virtual const int* get_dof_ids(mesh::TetMesh* m, mesh::Tet* tet) const = 0;
00064     virtual const bool* get_eliminated_dofs(mesh::TetMesh* m,
00065                                             mesh::Tet* tet,
00066                                             unsigned int sym_plane_config) const = 0;
00067 protected: // protected methods
00070     void get_Ae(mesh::Tet* tet, int order, const mesh::Tet::GradientMatrix& J, colarray::Matrix<double>& Ae);
00073     void get_Me(mesh::Tet* tet, int order, const mesh::Tet::GradientMatrix& J, colarray::Matrix<double>& Me);
00074 
00076     void constructGQ(int order);
00077     int            _intOrder;
00078     vector<double> _w;
00079     vector<double> _xi, _eta, _zeta;
00080 
00082     Matrix<double> _mu;
00083     Matrix<double> _eps;
00084 
00086     void evalElementFunctions(Matrix<double>& x, Matrix<double>& iM, 
00087                               Matrix<double>& Ned, Matrix<double>& Nod);
00088 
00090     Matrix<double> _Cur;
00091     Matrix<double> _Ned;
00092 
00095     void get_AeMe(mesh::Tet* tet, int order, colarray::Matrix<double>& Ae, colarray::Matrix<double>& Me);
00096 
00099     mesh::Vector3 eval_cartesian(mesh::Tet* tet, int order, const colarray::Vector<double>& dof, const mesh::Tet::GradientMatrix& J, const mesh::Vector3& x);
00100 
00103     mesh::Vector3 eval_curl_cartesian(mesh::Tet* tet, int order, const colarray::Vector<double>& dof, const mesh::Tet::GradientMatrix& J, const mesh::Vector3& x);
00104 
00107     mesh::Vector3 eval_simplex(mesh::Tet* tet, int order, const colarray::Vector<double>& dof, const mesh::Tet::GradientMatrix& J, const mesh::Vector4& x);
00111     mesh::Vector3 eval_curl_simplex(mesh::Tet* tet, int order, const colarray::Vector<double>& dof, const mesh::Tet::GradientMatrix& J, const mesh::Vector4& x);
00112 public:
00115     double surface_integral_curl(mesh::Tet* tet, mesh::Face* face, const colarray::Vector<double>& dofs);
00116 };
00117 
00118 #endif

Generated on Fri Oct 26 13:35:12 2007 for FEMAXX (Finite Element Maxwell Eigensolver) by  doxygen 1.4.7