VarFormPoisson Class Reference

#include <varformpoisson.h>

Inheritance diagram for VarFormPoisson:
Inheritance graph
[legend]
Collaboration diagram for VarFormPoisson:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 VarFormPoisson (FEOrthoMesh &fe)
 ~VarFormPoisson ()
virtual void assembly_local_system (OrthoMesh::Cell_It &cell, Matrix &ML)
virtual void assembly_local_rhs (OrthoMesh::Cell_It &cell, VecDouble &BL)
virtual const Matrixget_local_sparsity ()
virtual void mesh_change (OrthoMesh &mesh)
virtual void neumann_condition (VecDouble &Bl, OrthoMesh::Cell_It cell, OrthoMesh::Face_It face, unsigned face_dir, Function3D &fN)
void setParameters (VecDouble &K)

Private Attributes

Matrix _MSparsity
MappingOrthoMesh _map
MappedFEValues _feValues
std::vector< MappedFEFaceValues * > _vFaceValues
VecDoublepvK

Detailed Description

Definition at line 13 of file varformpoisson.h.


Constructor & Destructor Documentation

VarFormPoisson::VarFormPoisson ( FEOrthoMesh fe  ) 

Definition at line 4 of file varformpoisson.cpp.

00005    :VarFormOrthoMesh(fe),_feValues(fe)
00006  {
00007    assert(_fe.n_fields()==1);
00008    _MSparsity.reinit(_fe.n_local_dofs(),_fe.n_local_dofs());
00009    double c=1.0;
00010    _MSparsity.fill(c);
00011    dealii::QGauss<3> quad(2);
00012    dealii::QGauss<2> quad2d(2);
00013    
00014    _feValues.setPoints(quad);
00015    pvK=NULL;
00016    _vFaceValues.resize(UnitCube::n_faces());
00017 
00018    //set the fevaluesface for each face
00019    for (unsigned face_dir=0;face_dir<_vFaceValues.size();face_dir++)
00020    {
00021      _vFaceValues[face_dir]=new MappedFEFaceValues(fe,face_dir);
00022      _vFaceValues[face_dir]->setQuadrature(quad2d);
00023    }
00024  }

VarFormPoisson::~VarFormPoisson (  )  [inline]

Definition at line 25 of file varformpoisson.h.

00025 {}


Member Function Documentation

void VarFormPoisson::assembly_local_rhs ( OrthoMesh::Cell_It cell,
VecDouble BL 
) [virtual]

Assembly the local system right hand side not including the Neumman Condition. That is computed in another method neumann_condition()

Parameters:
cell The cell in the mesh to build the local system
BL To contain the right hand side of the local system The method assumes that the vector has already the correct size. You can get the local system rank using the dimension of the matrix returned by get_local_sparsity() method.

Implements VarFormOrthoMesh.

Definition at line 47 of file varformpoisson.cpp.

00048 {
00049   BL=0.0;
00050 }

void VarFormPoisson::assembly_local_system ( OrthoMesh::Cell_It cell,
Matrix ML 
) [virtual]

Assembly the local system

Parameters:
cell The cell in the mesh to build the local system
ML Matrix to contain the local system. The method assumes that the local matrix has already the correct sizes. You can get the local system rank using the dimension of the matrix returned by get_local_sparsity() method
Returns:
The local system at parameter ML

Implements VarFormOrthoMesh.

Definition at line 26 of file varformpoisson.cpp.

00027 {
00028   assert(pvK);
00029   assert(pvK->size() == getMesh()->numCells());
00030   unsigned n_local_dofs = _fe.n_local_dofs();
00031   unsigned n_pts = _feValues.n_pts();
00032   ML=0.0;
00033   double K = (*pvK)(cell->index());
00034   for (unsigned i=0;i<n_local_dofs;i++)
00035   {
00036     for (unsigned j=0;j<n_local_dofs;j++)
00037     {
00038       for (unsigned ptId=0;ptId<n_pts;ptId++)
00039       {
00040         ML(i,j)+=_feValues.function_grad(ptId,i)*_feValues.function_grad(ptId,j)*_feValues.JxW(ptId);
00041       }
00042     }
00043   }
00044   ML*=K;
00045 }

virtual const Matrix& VarFormPoisson::get_local_sparsity (  )  [inline, virtual]

Get the sparsity of the local system. The sparsity is represented by a full matrix with the same size as the local system matrix. A zero entry on position (i,j) means that this position is always zero in the local system. 1 means that such entry is always different of zero.

Returns:
The sparsity pattern of the local system

Implements VarFormOrthoMesh.

Definition at line 28 of file varformpoisson.h.

00028 {return _MSparsity;}

void VarFormPoisson::mesh_change ( OrthoMesh mesh  )  [virtual]

Implements VarFormOrthoMesh.

Definition at line 54 of file varformpoisson.cpp.

00055 {
00056   _map.reinit(mesh);
00057   _feValues.compute_values();
00058   _feValues.compute_grads(_map);
00059   _feValues.compute_JxW(_map);
00060 
00061   for (unsigned i=0;i<UnitCube::n_faces();i++)
00062   {
00063     _vFaceValues[i]->compute_values();
00064     //_vFaceValues[i]->compute_grads(_map);
00065     _vFaceValues[i]->compute_JxW(_map);
00066   }
00067   
00068 }

void VarFormPoisson::neumann_condition ( VecDouble Bl,
OrthoMesh::Cell_It  cell,
OrthoMesh::Face_It  face,
unsigned  face_dir,
Function3D fN 
) [virtual]

Set the RHS of the local system with the Neumman condition

Parameters:
Bl To contain the RHS with neumann condition This vector should be initialized with the right size
cell Cell of the local system
face Face where to get the neumman condition. The face must belong to the cell passed on the previous parameter.
face_id The direction of the face in respect to the cell Ex: LEFT_FACE, RIGHT_FACE, BOTTOM_FACE .....
fN The function that specifies the Neumman condition. The meaning of this function is specific to each variational formulation. Programmers Note: Always write a piece of documentation containing the meaning (what the value is) of the neumman function. For example FN(x) = Vel(x)*n(x), where n is the outward normal or n is the inward normal and so on.

Implements VarFormOrthoMesh.

Definition at line 72 of file varformpoisson.cpp.

00073 {
00074   //Compute the points to evaluate fN
00075   MappedFEFaceValues &fV = *( _vFaceValues[face_dir]);
00076   _map.reinit(cell);
00077   fV.compute_points(_map);
00078 
00079   unsigned n_local_dofs = _fe.face_to_cell_local_map(face_dir).size();
00080   unsigned n_pts = fV.n_pts();
00081   Bl=0.0;
00082   for (unsigned i=0;i<n_local_dofs;i++)
00083     for (unsigned ptId=0;ptId<n_pts;ptId++)
00084     {
00085       Bl(i)+=fV.function_value(ptId,i)*
00086         fN(fV.mapped_pts(ptId))*fV.JxW(ptId);
00087     }
00088 
00089 }

void VarFormPoisson::setParameters ( VecDouble K  )  [inline]

Definition at line 32 of file varformpoisson.h.

00032 {pvK=&K;}


Member Data Documentation

Definition at line 18 of file varformpoisson.h.

Definition at line 17 of file varformpoisson.h.

Definition at line 16 of file varformpoisson.h.

Definition at line 19 of file varformpoisson.h.

Definition at line 20 of file varformpoisson.h.


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