AssemblyOrthoMeshSystem< EqVarType, Matriz, Vetor > Class Template Reference

#include <assemblyorthomeshsystem.h>

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

List of all members.

Public Member Functions

 AssemblyOrthoMeshSystem ()
 ~AssemblyOrthoMeshSystem ()
void buildGlobalMatrix (Matriz &GM, OrthoMesh &mesh, EqVarType &eqvar, const bool checkZero=false, double tolerance=0.0)
template<class SparsePattern >
void buildSparsityPattern (SparsePattern &pattern, OrthoMesh &mesh, EqVarType &eqvar)
void buildGlobalRhs (Vetor &vRhs, OrthoMesh &mesh, EqVarType &eqvar)
void addNeummanCondition (const OrthoMesh &mesh, Vetor &vRHS, EqVarType &eqvar, Function3D &fN)
void setBoundaryValue (const OrthoMesh &mesh, EqVarType &eqvar, Function3D &f, MapIntDouble &map_boundary)

Protected Member Functions

void distribute_local_to_global (Matriz &GM, FullMatrix< double > &lM, const std::vector< unsigned > &local_dof_indices, const bool checkZero, double tolerance=0.0)

Private Attributes

FullMatrix< double > m_Ml
Vector< double > m_Bl

Detailed Description

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

Definition at line 12 of file assemblyorthomeshsystem.h.


Constructor & Destructor Documentation

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

Definition at line 57 of file assemblyorthomeshsystem.h.

00057 {}

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

Definition at line 58 of file assemblyorthomeshsystem.h.

00058 {}


Member Function Documentation

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystem< EqVarType, Matriz, Vetor >::addNeummanCondition ( const OrthoMesh mesh,
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 192 of file assemblyorthomeshsystem.h.

00193   {
00194     Point3D pt;
00195     const unsigned int   dofs_per_cell = eqvar.n_dofs_per_cell();//Number of dofs per cell
00196     const unsigned int n_face_quadrature_pts = eqvar.n_face_quadrature_points();//Number of QPoints per face
00197     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00198     endc = mesh.end_cell();
00199     
00200     for (; cell!=endc; cell++)//For all cells
00201     {
00202       if (cell->at_boundary()) //Verify if cell is at boundary
00203       {
00204         eqvar.setCell(cell);//If it is, set the cell and find faces at boundary
00205 
00206         //Iterate in all faces of the cell finding the cells at boundary
00207         for (unsigned face_no=0;face_no<OrthoMesh::FACES_PER_CELL;face_no++)
00208         {
00209           OrthoMesh::Face_It face = cell->face(static_cast<FaceDirection3D>(face_no));
00210           if (face->at_boundary()) //found face at boundary  
00211           {
00212             if (fN.isInDomain(face->barycenter()))
00213             {
00214               eqvar.setFace(face); //Set the face;
00215               //For all local dofs of the cell
00216               for (Index local_dof=0;local_dof<dofs_per_cell;local_dof++)
00217               {
00218                 //get what is the dof number of the face with respect to its local number
00219                 //in terms of the cell
00220                 Index cell_dof = eqvar.local_to_global_dofs_map()[local_dof]; 
00221                 double acum=0;
00222                 for (unsigned fQPoint=0;fQPoint < n_face_quadrature_pts;fQPoint++)
00223                 {
00224                   vRHS(cell_dof)+=eqvar.localNeummanCondition(face_no,local_dof,fQPoint,fN);
00225                   acum+=eqvar.localNeummanCondition(face_no,local_dof,fQPoint,fN);
00226                 }
00227                 printf("Acum: %g\n",acum);
00228               }
00229             }
00230           }
00231         }
00232       }
00233     }
00234   }

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystem< EqVarType, Matriz, Vetor >::buildGlobalMatrix ( Matriz &  GM,
OrthoMesh mesh,
EqVarType &  eqvar,
const bool  checkZero = false,
double  tolerance = 0.0 
) [inline]

Definition at line 60 of file assemblyorthomeshsystem.h.

00061   {
00062 
00063 
00064     unsigned i,j;
00065     const unsigned int   dofs_per_cell = eqvar.n_dofs_per_cell();
00066     const unsigned n_quadrature_points = eqvar.n_quadrature_points();
00067 
00068     //Inicializar matriz global e local.
00069     GM = 0;
00070     m_Ml.reinit(dofs_per_cell,dofs_per_cell);
00071    
00072     
00073 
00074     //Obter referencia de alguns campos de EqVarT
00075  
00076 
00077     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00078     endc = mesh.end_cell();
00079 
00080     for (; cell!=endc; cell++)
00081     {
00082       eqvar.setCell(cell);
00083       eqvar.changedCellEvent();
00084       //Build local matrix.
00085       m_Ml=0.0;
00086       for (Index qPoint=0; qPoint<n_quadrature_points; qPoint++)
00087       {
00088         for (i=0; i<dofs_per_cell; ++i)
00089         {
00090           for (j=0; j<dofs_per_cell; ++j)
00091             m_Ml(i,j) += eqvar.localEqVar(i,j,qPoint);
00092         }
00093       }
00094       //Distribua os nos na matriz global
00095       this->distribute_local_to_global(GM,m_Ml,eqvar.local_to_global_dofs_map(),checkZero,tolerance);
00096     }
00097     return;
00098   }

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystem< EqVarType, Matriz, Vetor >::buildGlobalRhs ( Vetor &  vRhs,
OrthoMesh mesh,
EqVarType &  eqvar 
) [inline]

Definition at line 140 of file assemblyorthomeshsystem.h.

00141     {
00142       unsigned i;
00143 
00144       const unsigned int   n_dofs_per_cell = eqvar.n_dofs_per_cell();
00145       unsigned n_quadrature_points = eqvar.n_quadrature_points();
00146 
00147 
00148       assert(vRhs.size() == eqvar.n_dofs());
00149       vRhs = 0;
00150 
00151 
00152 
00153       m_Bl.reinit(n_dofs_per_cell);
00154 
00155 
00156 
00157       OrthoMesh::Cell_It cell = mesh.begin_cell(),
00158       endc = mesh.end_cell();
00159 
00160       for (; cell!=endc; cell++)
00161       {
00162         eqvar.setCell(cell);
00163         m_Bl = 0;
00164         for (unsigned qPoint=0; qPoint<n_quadrature_points; qPoint++)
00165         {
00166           for (i=0; i<n_dofs_per_cell; ++i)
00167           {
00168             m_Bl(i) += eqvar.localRHS(i,qPoint);
00169           }
00170         }
00171         //Distribua os nos na matriz global
00172         for (unsigned int i=0; i<n_dofs_per_cell; ++i)
00173         {
00174           vRhs(eqvar.local_to_global_dofs_map()[i]) += m_Bl(i);
00175         }
00176       }
00177     }

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

Definition at line 102 of file assemblyorthomeshsystem.h.

00103   {
00104     unsigned i,j;
00105 
00106     //Initialize local matrix.
00107     m_Ml.reinit(eqvar.n_dofs_per_cell(),eqvar.n_dofs_per_cell());
00108    
00109     OrthoMesh::Cell_It cell = mesh.begin_cell(),
00110     endc = mesh.end_cell();
00111 
00112       //Build local matrix.
00113       m_Ml=0.0;
00114     unsigned n_dofs_per_cell=eqvar.n_dofs_per_cell();
00115 
00116     for (i=0; i<n_dofs_per_cell; ++i)
00117     {
00118       for (j=0; j<n_dofs_per_cell; ++j)
00119         m_Ml(i,j) += (eqvar.nullProduct(i,j)) ? 0 : 1;
00120     }
00121 
00122 
00123     for (; cell!=endc; cell++)
00124     {
00125       eqvar.setCell(cell);
00126       for (unsigned int i=0; i<n_dofs_per_cell; ++i)
00127         for (unsigned int j=0; j<n_dofs_per_cell; ++j)
00128         {
00129           if (m_Ml(i,j) == 0.0)
00130           {
00131             //printf("Not adding position (%d,%d) for value %g\n",local_dof_indices[i],local_dof_indices[j],m_Ml(i,j));
00132             continue;
00133           }
00134           else
00135             pattern.add (eqvar.local_to_global_dofs_map()[i],eqvar.local_to_global_dofs_map()[j]);
00136         }
00137     }
00138   }

template<class EqVarType, class Matriz = SparseMatrix<double>, class Vetor = VecDouble>
void AssemblyOrthoMeshSystem< EqVarType, Matriz, Vetor >::distribute_local_to_global ( Matriz &  GM,
FullMatrix< double > &  lM,
const std::vector< unsigned > &  local_dof_indices,
const bool  checkZero,
double  tolerance = 0.0 
) [inline, protected]

Definition at line 19 of file assemblyorthomeshsystem.h.

00020   {
00021     const unsigned dofs_per_cell = local_dof_indices.size();
00022     if (checkZero)
00023     {
00024       //Distribute the global matrix`s nodes. 
00025       for (unsigned int i=0; i<dofs_per_cell; ++i)
00026         for (unsigned int j=0; j<dofs_per_cell; ++j)
00027         {
00028           if (fabs(lM(i,j))<=tolerance)
00029             continue;
00030           else
00031           {
00032 #ifdef DEBUG
00033             if (GM.get_sparsity_pattern().exists(local_dof_indices[i],local_dof_indices[j]))
00034             {
00035               GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j));
00036             }
00037             else
00038               printf("Trying to put value %g in invalid entry (%d,%d)\n",lM(i,j),local_dof_indices[i],local_dof_indices[j]);
00039 #else
00040               GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j));
00041 #endif  
00042           }
00043         }
00044     }
00045     else
00046     {
00047       //Distribute the local matrix`s nodes. 
00048       for (unsigned int i=0; i<dofs_per_cell; ++i)
00049         for (unsigned int j=0; j<dofs_per_cell; ++j)
00050             GM.add (local_dof_indices[i],local_dof_indices[j],lM(i,j));
00051     }
00052 
00053   }

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

Definition at line 237 of file assemblyorthomeshsystem.h.

00238   {
00239     Point3D pt;
00240       const unsigned int   dofs_per_face = eqvar.n_dofs_per_face();
00241     std::vector<unsigned> face_dofs(dofs_per_face);
00242   OrthoMesh::Cell_It cell = mesh.begin_cell(),
00243         endc = mesh.end_cell();
00244     
00245     for (; cell!=endc; cell++)//For all cells
00246     {
00247       if (cell->at_boundary()) //Verify if cell is at boundary
00248       {
00249         eqvar.setCell(cell);//If it is, set the cell and find faces at boundary
00250 
00251         for (unsigned face_no=0;face_no<OrthoMesh::FACES_PER_CELL;face_no++)
00252         {
00253           
00254           if (cell->face(static_cast<FaceDirection3D>(face_no))->at_boundary())
00255           {
00256             if (f.isInDomain(cell->face(static_cast<FaceDirection3D>(face_no))->barycenter()))
00257             //For all local dofs at face
00258             for (Index local_face_dof=0;local_face_dof<dofs_per_face;local_face_dof++)
00259             {
00260               //get what is the dof number of the face with respect to its local number
00261               //in terms of the cell
00262               Index cell_dof = eqvar.local_face_to_local_cell_map(face_no,local_face_dof); 
00263               pt=eqvar.get_global_dof_position(cell_dof);
00264               map_boundary[eqvar.local_to_global_dofs_map()[cell_dof]]=f(pt,eqvar.local_dof_component(cell_dof));
00265               
00266             }
00267           }
00268         }
00269       }
00270     }
00271   }


Member Data Documentation

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

Definition at line 16 of file assemblyorthomeshsystem.h.

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

Definition at line 15 of file assemblyorthomeshsystem.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