PoissonFEM Class Reference

#include <poissonfem.h>

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

List of all members.

Public Member Functions

 PoissonFEM (OrthoMesh &mesh, VecDouble &cK, Function3D &fD, Function3D &fN, LinearSolver &solver, unsigned debugLevel=0)
 ~PoissonFEM ()
virtual void iterate (TransportBase &trans)
virtual void printOutput ()
virtual void getVelocitiesAtFaces (Matrix &M)
virtual const VecDoublegetPressureAtCells ()

Private Attributes

OrthoMeshm_mesh
VecDoublem_cK
FEDealWrapperQ1 ff
VarFormPoisson m_eqvar
SparsityPattern m_sparsity
SparseMatrix< double > m_M
VecDouble m_vRHS
VecDouble m_sol
LinearSolverm_solver
unsigned m_debugLevel
MapIntDouble m_dirichlet_bc
Function3Dm_fN

Detailed Description

PoissonFEM

Definition at line 15 of file poissonfem.h.


Constructor & Destructor Documentation

PoissonFEM::PoissonFEM ( OrthoMesh mesh,
VecDouble cK,
Function3D fD,
Function3D fN,
LinearSolver solver,
unsigned  debugLevel = 0 
)

Definition at line 11 of file poissonfem.cpp.

00012   :m_mesh(mesh), m_cK(cK),m_eqvar(ff),m_solver(solver),m_debugLevel(debugLevel), m_fN(fN)
00013 {
00014   m_eqvar.reinit(mesh);
00015   m_eqvar.setParameters(cK);
00016   AssemblyOrthoMeshSystemII<VarFormPoisson> ass;
00017   std::cout << "Resp 2: " << ass.max_nz_per_row(m_eqvar) << std::endl;
00018   unsigned rank = ass.global_system_rank(m_eqvar);
00019   std::cout << "Resp 1: " << rank << std::endl;
00020 
00021   //Build the sparsity patter and initialize the sparse matrix
00022   ass.buildSparsityPattern(m_sparsity,m_eqvar);
00023   m_M.reinit(m_sparsity);
00024   m_vRHS.reinit(rank);
00025   m_sol.reinit(rank);
00026   //Print the sparsity
00027   //std::ofstream file("sparse.dat");
00028   //m_sparsity.print_gnuplot(file);
00029 
00030   ass.setBoundaryValue(m_mesh,m_eqvar,fD,m_dirichlet_bc);
00031 
00032   /*
00033   AssemblyOrthoMeshSystem<PoissonEqVar> assembler;
00034   m_vRHS.reinit(m_eqvar.n_dofs());
00035   m_sol.reinit(m_eqvar.n_dofs());
00036   m_sparsity.reinit(m_eqvar.n_dofs(),m_eqvar.n_dofs(),m_eqvar.n_max_dof_coupling());
00037   assembler.buildSparsityPattern(m_sparsity,m_mesh,m_eqvar);
00038   m_sparsity.compress();
00039   m_M.reinit(m_sparsity);
00040   assembler.setBoundaryValue(m_mesh,m_eqvar,fD,dirichlet_bc_map);
00041   m_eqvar.setSource(*(new FConst3D(0)));
00042 
00043   */
00044 }

PoissonFEM::~PoissonFEM (  )  [inline]

Definition at line 34 of file poissonfem.h.

00034 {}


Member Function Documentation

virtual const VecDouble& PoissonFEM::getPressureAtCells (  )  [inline, virtual]

Reimplemented from DynamicBase.

Definition at line 38 of file poissonfem.h.

00038 {return m_sol;}

void PoissonFEM::getVelocitiesAtFaces ( Matrix M  )  [virtual]

Reimplemented from DynamicBase.

Definition at line 112 of file poissonfem.cpp.

00113 {
00114   throw new Exception("PoissonFEM::getVelocitiesAtFaces Not implemented yet\n");
00115 }

void PoissonFEM::iterate ( TransportBase trans  )  [virtual]

Implements DynamicBase.

Definition at line 47 of file poissonfem.cpp.

00048 {
00049   
00050   printf("Building system\n");
00051   AssemblyOrthoMeshSystemII<VarFormPoisson> assembler;
00052   assembler.buildGlobalMatrix(m_M, m_eqvar);
00053   assembler.buildGlobalRHS(m_vRHS,m_eqvar);
00054 
00055   assembler.addNeumannCondition(m_vRHS,m_eqvar,m_fN);
00056   //Apply boundary conditions
00057   m_vRHS.print(std::cout);
00058   assembler.applyBoundaryValue(m_dirichlet_bc,m_vRHS);
00059   assembler.applyBoundaryValue(m_dirichlet_bc,m_M);
00060 
00061   //solve the system
00062   printf("Solving\n");
00063   m_solver.solve(m_M,m_sol,m_vRHS);
00064   printf("done");
00065   
00066   //std::ofstream file("matrix.dat");
00067   //m_M.print_formatted(file,10,true,0," 0 ");
00068   //m_vRHS.print(std::cout);
00069   /*
00070   m_eqvar.setK(m_cK);
00071   assembler.buildGlobalRhs(m_vRHS,m_mesh,m_eqvar);
00072   assembler.applyBoundaryValue(dirichlet_bc_map,m_M);
00073   assembler.applyBoundaryValue(dirichlet_bc_map,m_vRHS);
00074   assembler.addNeummanCondition(m_mesh,m_vRHS,m_eqvar,m_fN);
00075   m_sol=m_vRHS;
00076   m_sol.print();
00077   solveA.factorize(m_M);
00078   solveA.solve(m_sol);
00079   */
00080 }

void PoissonFEM::printOutput (  )  [virtual]

Implements DynamicBase.

Definition at line 84 of file poissonfem.cpp.

00085 {
00086   AssemblyOrthoMeshSystemII<VarFormPoisson> assembler;
00087   VecDouble dofSol(m_sol.size());
00088   assembler.getDOFSolAtVertices(dofSol,m_sol,m_eqvar);
00089   HDF5OrthoWriter &hdf5 = HDF5OrthoWriter::getHDF5OrthoWriter(); 
00090   hdf5.writeScalarField(dofSol,"P");
00091 
00092   //dofSol-=m_sol;
00093 
00094   /*
00095   if (m_debugLevel > 0)
00096   {
00097     HDF5OrthoWriter &hdf5 = HDF5OrthoWriter::getHDF5OrthoWriter(); 
00098     hdf5.writeScalarField(m_sol,"P");
00099 
00100      GnuPlotAnim &plot = GnuPlotAnim::getGnuPlotAnim();
00101      
00102      //char str[100];
00103      //printf("Mass: %g\n",m_mesh.getIntegralAtCells(m_cValues));
00104      plot.set("yrange","[0:1]");
00105      plot.plotVerticeValuesAtFixedZ(m_mesh,m_sol,"Sc","with lines",5);
00106      plot.NextScene();
00107   }
00108   */
00109 }


Member Data Documentation

Definition at line 20 of file poissonfem.h.

Definition at line 19 of file poissonfem.h.

unsigned PoissonFEM::m_debugLevel [private]

Definition at line 26 of file poissonfem.h.

Definition at line 27 of file poissonfem.h.

Definition at line 21 of file poissonfem.h.

Definition at line 28 of file poissonfem.h.

SparseMatrix<double> PoissonFEM::m_M [private]

Definition at line 23 of file poissonfem.h.

Definition at line 18 of file poissonfem.h.

Definition at line 24 of file poissonfem.h.

Definition at line 25 of file poissonfem.h.

SparsityPattern PoissonFEM::m_sparsity [private]

Definition at line 22 of file poissonfem.h.

Definition at line 24 of file poissonfem.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:28 2012 for CO2INJECTION by  doxygen 1.6.3