AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor > Class Template Reference

#include <assemblyorthomeshsystemII.h>

Inheritance diagram for AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >:
Inheritance graph
[legend]
Collaboration diagram for AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 AssemblyOrthoMeshSystemII ()
 ~AssemblyOrthoMeshSystemII ()
void buildGlobalMatrix (Matriz &GM, EqVarType &eqvar)
void buildGlobalRHS (Vetor &vRHS, EqVarType &eqvar)
Index global_system_rank (EqVarType &eqvar)
Index max_nz_per_row (EqVarType &eqvar)
template<class SparsePattern >
void buildSparsityPattern (SparsePattern &pattern, EqVarType &eqvar)
void addNeumannCondition (Vetor &vRHS, EqVarType &eqvar, Function3D &fN)
void setBoundaryValue (const OrthoMesh &mesh, EqVarType &eqvar, Function3D &f, MapIntDouble &map_boundary)
void getDOFSolAtVertices (VecDouble &vSol, const VecDouble &dofSol, EqVarType &eqvar, unsigned component=0)

Private Attributes

Matrix m_Ml
VecDouble m_Bl

Detailed Description

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
class AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >

Definition at line 14 of file assemblyorthomeshsystemII.h.


Constructor & Destructor Documentation

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::AssemblyOrthoMeshSystemII (  )  [inline]

Definition at line 18 of file assemblyorthomeshsystemII.h.

00022 :

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::~AssemblyOrthoMeshSystemII (  )  [inline]

Definition at line 19 of file assemblyorthomeshsystemII.h.

00022 :


Member Function Documentation

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::addNeumannCondition ( Vetor &  vRHS,
EqVarType &  eqvar,
Function3D fN 
) [inline]

Adds neumman condition to the vector. The condition is inserted using the function fN defined in the boundaries of the domain. For each face at boundary this method call addNeummanCondition method of the EqVarType class passing the local dof of the face, que quadrature point and the function fN. Is the job of the EqVarType to return the correct neumman condition for each test function corresponding to the local dof face and quadrature points

Parameters:
mesh Mesh
vRhs The right hand side vector
eqvar Variational equation class
fN Function defining the Neumman condition
Returns:

Definition at line 236 of file assemblyorthomeshsystemII.h.

00244   {
00245     assert(eqvar.getMesh());
00246     OrthoMesh &mesh=*(eqvar.getMesh());
00247     
00248     Point3D pt;
00249     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00250     endc = mesh.end_cell();
00251 
00252     VecDouble Bl(eqvar.getFE().n_dofs_per_face());
00253     
00254     for (; cell!=endc; cell++)//For all cells
00255     {
00256       if (cell->at_boundary()) //Verify if cell is at boundary
00257       {
00258         eqvar.reinit(cell);//If it is, set the cell and find faces at boundary
00259 
00260         //Iterate in all faces of the cell finding the cells at boundary
00261         for (unsigned face_dir=0;face_dir<OrthoMesh::FACES_PER_CELL;face_dir++)
00262         {
00263           OrthoMesh::Face_It face = cell->face(static_cast<FaceDirection3D>(face_dir));
00264           if (face->at_boundary()) //found face at boundary  
00265           {
00266             if (fN.isInDomain(face->barycenter()))
00267             {
00268               const VecIndex& global_map = eqvar.getFE().get_global_map(cell);
00269               const VecIndex& face_local_map = eqvar.getFE().face_to_cell_local_map(face_dir);
00270 
00271               //get the neuwmann condition
00272               eqvar.neumann_condition(Bl,cell,face,face_dir,fN);
00273 
00274 
00275               //For all local dofs of the face
00276               for (unsigned face_dof=0;face_dof<face_local_map.size();face_dof++)
00277               {
00278                 //get what is the dof number of the face with respect to its local number
00279                 //in terms of the cell
00280                 unsigned local_cell_dof=face_local_map[face_dof];
00281                 vRHS(global_map[local_cell_dof])+=Bl(face_dof);
00282 

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::buildGlobalMatrix ( Matriz &  GM,
EqVarType &  eqvar 
) [inline]

Definition at line 21 of file assemblyorthomeshsystemII.h.

00022            :
00023 
00024 
00025   public:
00026     AssemblyOrthoMeshSystemII(){}
00027     ~AssemblyOrthoMeshSystemII(){}
00028   
00029   void buildGlobalMatrix(Matriz& GM, EqVarType& eqvar)
00030   {
00031     GM=0;
00032     assert(eqvar.getMesh());
00033     OrthoMesh &mesh=*(eqvar.getMesh());
00034 
00035 
00036     //Get the sparsity pattern.
00037     const Matrix &MSparsity = eqvar.get_local_sparsity();
00038 
00039     assert(MSparsity.m() == MSparsity.n());
00040     unsigned n_rank = MSparsity.m();
00041     //Initialize  global and local matrix.
00042     //The size of the local system is equal to the MSparsity
00043     //rank
00044     m_Ml.reinit(n_rank,n_rank);
00045     //For all cells
00046     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00047     endc = mesh.end_cell();
00048     for (; cell!=endc; cell++)
00049     {
00050       eqvar.reinit(cell);
00051       const VecIndex &local_global_map = eqvar.getFE().get_global_map(cell);
00052 
00053       eqvar.assembly_local_system(cell,m_Ml);
00054 
00055 
00056       //Distribute the local matrix into the global matrix
00057       for (unsigned i=0;i<n_rank;i++)
00058       {
00059         for (unsigned j=0;j<n_rank;j++)
        {

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::buildGlobalRHS ( Vetor &  vRHS,
EqVarType &  eqvar 
) [inline]

Definition at line 61 of file assemblyorthomeshsystemII.h.

00069     {
00070 
00071     assert(eqvar.getMesh());
00072     OrthoMesh &mesh=*(eqvar.getMesh());
00073     vRHS=0;
00074 
00075     //Get the rank of the local system
00076     unsigned n_rank = eqvar.get_local_sparsity().m();
00077 
00078 
00079 
00080     //Initialize  local RHS
00081     m_Bl.reinit(n_rank);
00082 
00083     //For all cells
00084     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00085     endc = mesh.end_cell();
00086     for (; cell!=endc; cell++)
00087     {
00088       eqvar.reinit(cell);
00089       const VecIndex &local_global_map = eqvar.getFE().get_global_map(cell);
00090 
00091       //Build local rhs
00092       eqvar.assembly_local_rhs(cell,m_Bl);
00093 
00094 
00095       //Distribute the local matrix into the global matrix
00096 
00097       for (unsigned i=0;i<n_rank;i++)

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
template<class SparsePattern >
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::buildSparsityPattern ( SparsePattern &  pattern,
EqVarType &  eqvar 
) [inline]

Definition at line 179 of file assemblyorthomeshsystemII.h.

00187   {
00188 
00189     assert(eqvar.getMesh());
00190     OrthoMesh &mesh=*(eqvar.getMesh());
00191 
00192     //Get the global system ranks
00193     Index global_rank = global_system_rank(eqvar);
00194 
00195     //Initialize the sparse pattern calculating the
00196     //max non zero per rows with the function max_nz_per_row()
00197     pattern.reinit(global_rank,
00198                    global_rank,
00199                    this->max_nz_per_row(eqvar));
00200     
00201     //Get the local sparsity
00202     const Matrix &localSparsity = eqvar.get_local_sparsity();
00203     assert(localSparsity.m() == localSparsity.n());
00204     unsigned n_rank=localSparsity.m();
00205 
00206     
00207     //Now a loop to build the sparsity pattern
00208     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00209     endc = mesh.end_cell();
00210     for (; cell!=endc; cell++)
00211     {
00212       eqvar.reinit(cell);
00213       const VecIndex &local_global_map = eqvar.getFE().get_global_map(cell);
00214       for (unsigned i=0;i<n_rank;i++)
00215         for (unsigned j=0;j<n_rank;j++)
00216         {
00217           if (localSparsity(i,j) != 0.0)
00218             pattern.add(local_global_map[i],local_global_map[j]);
00219         }
00220     }

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::getDOFSolAtVertices ( VecDouble vSol,
const VecDouble dofSol,
EqVarType &  eqvar,
unsigned  component = 0 
) [inline]

Definition at line 339 of file assemblyorthomeshsystemII.h.

00347   {
00348     //Set a list of points in the vertices
00349     
00350     MappedFEValues feValues(eqvar.getFE());;
00351     UnitCube refDomain;
00352     VecDouble weights(refDomain.n_vertices());
00353     weights=1;
00354     feValues.setPoints(refDomain.vertices(),weights);
00355     feValues.compute_values();
00356     
00357 
00358     //Initially put the sol values with INFINITY value;
00359     vSol=INFINITY;
00360 
00361     //For all cells of the eqvar's mesh
00362     assert(eqvar.getMesh());
00363     OrthoMesh &mesh = *(eqvar.getMesh());
00364     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00365     endc = mesh.end_cell();
00366     for (; cell!=endc; cell++)//For all cells
00367     {
00368       //get the global map;
00369       const VecIndex &global_map=  eqvar.getFE().get_global_map(cell);
00370 
00371       //for all vertices
00372       for (unsigned v=0;v<refDomain.n_vertices();v++)
00373       {
00374         //Get vertice on the mesh
00375         unsigned vIndex=cell->vertex_index(static_cast<VertexDirection3D>(v));
00376         //Check if the solution was alread set
00377         if (vSol(vIndex) == INFINITY)
00378         {
00379           //Ok the position for this vertice was not set
00380           //do it.
00381           vSol(vIndex)=0;
00382           for (unsigned fId=0;fId<global_map.size();fId++)
00383           {

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
Index AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::global_system_rank ( EqVarType &  eqvar  )  [inline]

Definition at line 100 of file assemblyorthomeshsystemII.h.

00108   {
00109     
00110     Index max=0;
00111     assert(eqvar.getMesh() != NULL);
00112     OrthoMesh &mesh = *(eqvar.getMesh());
00113     //For all cells
00114     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00115     endc = mesh.end_cell();
00116     for (; cell!=endc; cell++)
00117     {
00118       eqvar.reinit(cell);
00119       const VecIndex &local_global_map = eqvar.getFE().get_global_map(cell);
00120       for (unsigned i=0;i<local_global_map.size();i++)

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
Index AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::max_nz_per_row ( EqVarType &  eqvar  )  [inline]

Definition at line 124 of file assemblyorthomeshsystemII.h.

00132   {
00133     Point3D P1(0,0,0),P2(1,1,1);
00134     
00135     //Create a small mesh
00136     OrthoMesh mesh(P1,P2,3,3,3);
00137 
00138     //Save the mesh already assigned to the eqvar
00139     OrthoMesh *backMesh = eqvar.getMesh();
00140 
00141     //Now initialize the finite element with the new mesh.
00142     eqvar.reinit(mesh);
00143 
00144     //Get the global system ranks
00145     //for the small system of a mesh 3x3x3
00146     Index global_rank = global_system_rank(eqvar);
00147 
00148     //Initialize the sparse pattern as a full sparse pattern
00149     //(i.e The number of non zero columns per row is the
00150     //number of columns)
00151     SparsityPattern pattern;
00152     pattern.reinit(global_rank,global_rank,global_rank);
00153     
00154     //Get the local sparsity
00155     const Matrix &localSparsity = eqvar.get_local_sparsity();
00156     assert(localSparsity.m() == localSparsity.n());
00157     unsigned n_rank=localSparsity.m();
00158     //Now a loop to build the sparsity pattern
00159     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00160     endc = mesh.end_cell();
00161     for (; cell!=endc; cell++)
00162     {
00163       eqvar.reinit(cell);
00164       const VecIndex &local_global_map = eqvar.getFE().get_global_map(cell);
00165       for (unsigned i=0;i<n_rank;i++)
00166         for (unsigned j=0;j<n_rank;j++)
00167         {
00168           if (localSparsity(i,j) != 0.0)
00169             pattern.add(local_global_map[i],local_global_map[j]);
00170         }
00171     }
00172     //Ok, now revert the finite element with the previous mesh
00173     eqvar.reinit(*backMesh);
00174 
00175 
00176     //Finally compress the pattern

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::setBoundaryValue ( const OrthoMesh mesh,
EqVarType &  eqvar,
Function3D f,
MapIntDouble map_boundary 
) [inline]

Definition at line 285 of file assemblyorthomeshsystemII.h.

00293   {
00294     assert(&mesh == eqvar.getMesh());
00295     VecDouble pt(3);
00296     VecDouble ptRef(3);
00297     MappingOrthoMesh map;
00298     map.reinit(mesh);
00299     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00300     endc = mesh.end_cell();
00301     for (; cell!=endc; cell++)//For all cells
00302     {
00303       if (cell->at_boundary()) //Verify if cell is at boundary
00304       {
00305         eqvar.reinit(cell);//If it is, set the cell 
00306 
00307         //For all faces of this cell
00308         for (unsigned face_no=0;face_no<OrthoMesh::FACES_PER_CELL;face_no++)
00309         {
00310           //if face is at boundary
00311           if (cell->face(static_cast<FaceDirection3D>(face_no))->at_boundary())
00312           {
00313             //See if the barycenter of the face belongs to the domain
00314             //of the boundary condition function
00315             if (f.isInDomain(cell->face(static_cast<FaceDirection3D>(face_no))->barycenter()))
00316             {
00317               map.reinit(cell);
00318               const VecIndex &global_map=eqvar.getFE().get_global_map(cell);
00319               const VecIndex &face_dofs= eqvar.getFE().face_to_cell_local_map(face_no);
00320               //For all local dofs at face
00321               for (Index local_face_dof=0;local_face_dof<face_dofs.size();local_face_dof++)
00322               {
00323                 //get what is the dof number of the face with respect to its local number
00324                 //in terms of the cell
00325                 Index cell_dof = face_dofs[local_face_dof]; 
00326 
00327                 //Get the support point for this dof in ptRef
00328                 ptRef=eqvar.getFE().get_support_points()[cell_dof];
00329 
00330                 //Map to the point in the  original domain
00331                 map.T(ptRef,pt);
00332                 
00333                 //Set the value of the dof with the value of the
00334                 //boundary condition function evaluated in the
00335                 //support point of the deformed domain
00336                 map_boundary[global_map[cell_dof]]=f(pt,eqvar.getFE().local_dof_field(cell_dof));


Member Data Documentation

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
VecDouble AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::m_Bl [private]

Definition at line 11 of file assemblyorthomeshsystemII.h.

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
Matrix AssemblyOrthoMeshSystemII< EqVarType, Matriz, Vetor >::m_Ml [private]

Definition at line 10 of file assemblyorthomeshsystemII.h.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Sun Apr 8 23:12:56 2012 for CO2INJECTION by  doxygen 1.6.3