Methods to get values associated with mesh geometries from functions

Functions

void DealBase::setCentralValuesFromFunction (Function3D &f, VecDouble &cValues, unsigned component=0)
void DealBase::setCentralValuesInFunctionDomain (Function3D &f, MapIntDouble &values, unsigned component=0)
void DealBase::setVerticesValuesFromFunction (Function< DIM > &f, VecDouble &vValues, unsigned deg=0)
void DealBase::setFacesValuesFromFunction (Function3D &f, VecDouble &fValues, unsigned deg=0, bool bOnlyPrescribedBoundary=false)
void DealBase::setFacesValuesFromFunction (Function3D &f, Matrix &MValues, bool bOnlyPrescribedBoundary=false)
void DealBase::tagFacesInDomain (Function3D &f, VecTag &tags, char value)
void DealBase::tagFacesInDomain (Function3D &f, VecBool &tags)
void DealBase::setVerticeValue (Cell &cell, VertexDirection3D vertex, double dd, VecDouble &values)
void DealBase::setMaterialIdOfPrescribedBoundaryCells (GeneralFunctionInterface &f, unsigned material_id)
void DealBase::setMaterialIdFromFunctionDomain (GeneralFunctionInterface &f, unsigned material_id)
void DealBase::getCellsInFunctionDomain (Function3D &f, VecIndex &vec)

Function Documentation

void DealBase::getCellsInFunctionDomain ( Function3D f,
VecIndex vec 
) [inherited]

Get a list of the indices of all cells whose barycenter is inside the domain of the function f.

Parameters:
f function
vec vector to contain the indices (output)

Definition at line 938 of file dealbase.cpp.

00939 {
00940   vec.clear();
00941   Cell cell = getTriangulation().begin_active(),
00942        endc = getTriangulation().end();
00943   for (;cell!=endc;++cell)
00944   {
00945     Point3D barycenter = cell->barycenter(); 
00946     if (f.isInDomain(barycenter))
00947     {
00948       vec.push_back(cell->index());
00949     }
00950     
00951   }
00952 }

void DealBase::setCentralValuesFromFunction ( Function3D f,
VecDouble cValues,
unsigned  component = 0 
) [inherited]

Eval the function f at the central values in the cell, storing the solution in the vector cValues.

Parameters:
f Function to be evaluated.
cValues Vector containing the central values.
component The component of the method.

Definition at line 753 of file dealbase.cpp.

00754 {
00755   assert(cValues.size() == numCells());
00756          
00757   Cell cell = getTriangulation().begin_active(),
00758        endc = getTriangulation().end();
00759   Point<DIM> BC;     
00760   for (;cell!=endc;++cell)
00761   {
00762     BC = cell->barycenter();
00763     setCentralValue(cell,f(BC,component),cValues);
00764   }
00765 }

void DealBase::setCentralValuesInFunctionDomain ( Function3D f,
MapIntDouble values,
unsigned  component = 0 
) [inherited]

Set values of a cell in a Map structure containing only the indices of the cells in the domain of the function

Parameters:
f Function
values Map containing the value of f for each cell in the domain of f
component The component of the function f
Returns:

Definition at line 773 of file dealbase.cpp.

00774 {
00775          
00776   Cell cell = getTriangulation().begin_active(),
00777        endc = getTriangulation().end();
00778   Point<DIM> BC;     
00779   for (;cell!=endc;++cell)
00780   {
00781     BC = cell->barycenter();
00782     if (f.isInDomain(BC,component))
00783       values[cell->index()] = f(BC,component);
00784   }
00785 }

void DealBase::setFacesValuesFromFunction ( Function3D f,
Matrix MValues,
bool  bOnlyAtBoundary = false 
) [inherited]

Eval all the components of the function f at the faces midpoints storing the values in the Matrix fValues.

Parameters:
f Function to be evaluated.
fValues Matrix to contain the values. The Matrix must have the size [numFaces() x f.num_components()]
bOnlyAtBoundary If the function must be evaluated only in the boundary faces.
Returns:

Definition at line 875 of file dealbase.cpp.

00876 {
00877   assert(MValues.m() == numFaces());
00878   assert(MValues.n() == f.n_components());
00879   
00880   Face face    = getTriangulation().begin_active_face();
00881   Face endFace = getTriangulation().end_face();
00882   
00883   for (;face!=endFace;face++)
00884   {
00885     if (bOnlyAtBoundary && !face->at_boundary()) 
00886       continue;
00887 
00888     Point<3> P = face->vertex(BL_VERTEX);
00889     P+=face->vertex(UR_VERTEX);
00890     P[0]=P[0]/2.0;
00891     P[1]=P[1]/2.0;
00892     P[2]=P[2]/2.0;
00893     unsigned faceIndex = face->index();
00894     for (unsigned deg=0;deg<f.n_components();deg++)
00895     {
00896       MValues(faceIndex,deg) = f(P,deg);
00897     }
00898   }
00899 }

void DealBase::setFacesValuesFromFunction ( Function3D f,
VecDouble fValues,
unsigned  deg = 0,
bool  bOnlyPrescribedBoundary = false 
) [inherited]

Eval the function f at the faces midpoints storing the values in the vector fValues.

Parameters:
f Function to be evaluated.
fValues Vector to contain the values.
deg The number of the component of the function f to be stored in fValues.
bOnlyAtBoundary If this flag is true the function is evaluated only in the boundary faces and only in the faces where the method f.boundary_type(Point X) returns the value PRESCRIBED_BOUNDARY_CONDITION.
Returns:

Definition at line 797 of file dealbase.cpp.

00798 {
00799   assert(fValues.size() == numFaces());
00800   Face face    = getTriangulation().begin_active_face();
00801   Face endFace = getTriangulation().end_face();
00802   
00803   for (;face!=endFace;face++)
00804   {
00805     if (bOnlyPrescribedBoundary && !face->at_boundary()) 
00806       continue;
00807     
00808     Point<3> P = face->vertex(BL_VERTEX);
00809     P+=face->vertex(UR_VERTEX);
00810     P[0]=P[0]/2.0;
00811     P[1]=P[1]/2.0;
00812     P[2]=P[2]/2.0;
00813     if (!f.isInDomain(P,deg))
00814       continue;
00815     else
00816       setFaceValue(face,f(P,deg),fValues);
00817   }
00818 }

void DealBase::setMaterialIdFromFunctionDomain ( GeneralFunctionInterface f,
unsigned  material_id 
) [inherited]

Objetivo: Configurar as celulas cujo baricentro esta no dominio de f (Ver Function2D::isInDomain(Point)

Parametros: f = Funcao cujo metoto isInDomain sera utilizado. material_id = Identificador do material a ser configurado na celula.

Definition at line 190 of file dealbase.cpp.

00191 {
00192   Cell cell = getTriangulation().begin_active();
00193   Cell endc = getTriangulation().end();
00194 
00195   for (;cell!=endc;cell++)
00196   {
00197     if (f.isInDomain(cell->barycenter()))
00198     {
00199       cell->set_material_id(material_id);
00200       cout << "done \n";
00201     }
00202 
00203   }
00204 }

void DealBase::setMaterialIdOfPrescribedBoundaryCells ( GeneralFunctionInterface f,
unsigned  material_id 
) [inherited]

Objetivo: Configurar materialId de todas as celulas que possui pelo menos um vertice de fronteira prescrito por f.

Parametros: f = Funcao de contorno. material_id

Definition at line 216 of file dealbase.cpp.

00217 {
00218   Triangulation<DIM> &tria = getTriangulation();
00219   Cell cell = tria.begin_active();
00220   Cell endc = tria.end();
00221   
00222 
00223   for (;cell!=endc;cell++)
00224   {
00225     if (cell->at_boundary())
00226     {
00227       for (unsigned i=0;i<GeometryInfo<DIM>::vertices_per_cell;i++)
00228       {
00229         if (f.isInDomain(cell->vertex(i)))
00230         {
00231           cell->set_material_id(material_id);
00232           break;
00233         }
00234       }
00235     }
00236   }
00237 }

void DealBase::setVerticesValuesFromFunction ( Function< DIM > &  f,
VecDouble vValues,
unsigned  deg = 0 
) [inherited]

Definition at line 77 of file dealbase.cpp.

00078 {
00079   assert(vValues.size() == numVertices());
00080   Cell cell = getTriangulation().begin_active(),
00081        endc = getTriangulation().end();
00082   for (;cell!=endc;++cell)
00083   {
00084     for (unsigned i=0;i<GeometryInfo<DIM>::vertices_per_cell;i++)
00085     {
00086       setVerticeValue(cell,(VertexDirection3D) i,f.value(cell->vertex(i),deg),vValues);
00087     }
00088   }
00089 }

void DealBase::setVerticeValue ( Cell cell,
VertexDirection3D  vertex,
double  dd,
VecDouble values 
) [inherited]

Objetivo:

Parametros:

Retorno:

Definition at line 248 of file dealbase.cpp.

00249 {
00250   assert(values.size() == numVertices());
00251   values(cell->vertex_index(vertex)) = dd;   
00252 
00253 }

void DealBase::tagFacesInDomain ( Function3D f,
VecBool tags 
) [inherited]

Definition at line 847 of file dealbase.cpp.

00848 {
00849   assert(tags.size() == numFaces());
00850   Face face    = getTriangulation().begin_active_face();
00851   Face endFace = getTriangulation().end_face();
00852   
00853   for (;face!=endFace;face++)
00854   {
00855     Point<3> P = face->vertex(BL_VERTEX);
00856     P+=face->vertex(UR_VERTEX);
00857     P[0]=P[0]/2.0;
00858     P[1]=P[1]/2.0;
00859     P[2]=P[2]/2.0;
00860     if (f.isInDomain(P,0))
00861       tags[face->index()] = true;
00862   }
00863 }

void DealBase::tagFacesInDomain ( Function3D f,
VecTag tags,
char  value 
) [inherited]

Tag faces inside the domain of f using the value value.

Parameters:
f Function
tags Vector of tags
value Value to tag
Returns:

Definition at line 828 of file dealbase.cpp.

00829 {
00830   assert(tags.size() == numFaces());
00831   Face face    = getTriangulation().begin_active_face();
00832   Face endFace = getTriangulation().end_face();
00833   
00834   for (;face!=endFace;face++)
00835   {
00836     Point<3> P = face->vertex(BL_VERTEX);
00837     P+=face->vertex(UR_VERTEX);
00838     P[0]=P[0]/2.0;
00839     P[1]=P[1]/2.0;
00840     P[2]=P[2]/2.0;
00841     if (f.isInDomain(P,0))
00842       tags[face->index()] = value;
00843   }
00844 }

 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