hp3D::Identity< F > Class Template Referenceabstract

A function class to calculate element matrices for the mass matrix. More...

#include <bilinearForm.hh>

Inheritance diagram for hp3D::Identity< F >:
concepts::BilinearForm< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
 Identity (const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
 Constructor. More...
 
virtual void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
 
virtual void operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em) const=0
 Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em. More...
 
virtual void operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em, const ElementPair< typename Realtype< Real >::type > &ep) const
 Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em. More...
 
void operator() (const Hexahedron &elmX, const Hexahedron &elmY, concepts::ElementMatrix< F > &em) const
 
virtual ~Identity ()
 

Timing Interface

These functions are used to get timings from class internal computations.

The values are stored in a user defined concepts::InOutParameters structure in different arrays (see setTimings). These arrays can be grouped into a table for easier postprocessing with

table.addMap(concepts::ResultsTable::DOUBLE, "jacobian", output);
table.addMap(concepts::ResultsTable::DOUBLE, "whole_sumfact", output);
std::ofstream ofs("table.gnuplot");
ofs << std::setprecision(20);
concepts::Array< F > jacobian_
 Intermediate data for element matrix computation. More...
 
const HexahedronoldElm_
 Intermediate data for element matrix computation. More...
 
concepts::ElementFormulaContainer< F > frm_
 Element formula. More...
 
static concepts::InOutParameterstimings_
 Place to store timing values. More...
 
static uint timeCntr_
 Counter for timing table. More...
 
virtual Identityclone () const
 Intermediate data for element matrix computation. More...
 
static void setTimings (concepts::InOutParameters *timings)
 Sets the class to store the timing values in. More...
 
static bool timings ()
 Returns true if the class is able to do timings. More...
 
virtual std::ostream & info (std::ostream &os) const
 Intermediate data for element matrix computation. More...
 

Detailed Description

template<class F = Real>
class hp3D::Identity< F >

A function class to calculate element matrices for the mass matrix.

This class is equiped with an interface to get timings of internal computations if compiled accordingly (see bilinearForm.cc file), see setTimings() and timings().

Author
Philipp Frauenfelder, 2001
Examples
hpFEM3d-EV.cc.

Definition at line 84 of file bilinearForm.hh.

Constructor & Destructor Documentation

◆ Identity()

template<class F = Real>
hp3D::Identity< F >::Identity ( const concepts::ElementFormulaContainer< F >  frm = concepts::ElementFormulaContainer< F >())

Constructor.

◆ ~Identity()

template<class F = Real>
virtual hp3D::Identity< F >::~Identity ( )
virtual

Member Function Documentation

◆ clone() [1/2]

template<class F = Real>
virtual Identity* hp3D::Identity< F >::clone ( ) const
inlinevirtual

Intermediate data for element matrix computation.

Implements concepts::Cloneable.

Definition at line 129 of file bilinearForm.hh.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< Real , typename Realtype<Real >::type >::clone ( ) const
pure virtualinherited

Virtual constructor.

Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.

◆ info()

template<class F = Real>
virtual std::ostream& hp3D::Identity< F >::info ( std::ostream &  os) const
protectedvirtual

Intermediate data for element matrix computation.

Reimplemented from concepts::BilinearForm< Real >.

◆ operator()() [1/4]

template<class F = Real>
virtual void hp3D::Identity< F >::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< F > &  em 
) const
virtual

◆ operator()() [2/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  em 
) const
pure virtualinherited

Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em.

Postcondition
The returned matrix em has the correct size.
Parameters
elmXLeft element (test functions)
elmYRight element (trial functions)
emReturn element matrix

◆ operator()() [3/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  em,
const ElementPair< G > &  ep 
) const
inlinevirtualinherited

Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em.

If this method is not reimplemented in a derived class, the default behaviour is to call the application operator without ep.

Postcondition
The returned matrix em has the correct size.
Parameters
elmXLeft element
elmYRight element
emReturn element matrix
epElement pair holding more information on the pair elmX and elmY

Definition at line 57 of file bilinearForm.hh.

◆ operator()() [4/4]

template<class F = Real>
void hp3D::Identity< F >::operator() ( const Hexahedron elmX,
const Hexahedron elmY,
concepts::ElementMatrix< F > &  em 
) const

◆ setTimings()

template<class F = Real>
static void hp3D::Identity< F >::setTimings ( concepts::InOutParameters timings)
static

Sets the class to store the timing values in.

Additionally, the timeCntr_ is reset to 0. This counter is used to fill in the values into the arrays listed below in subsequent calls. The following timings are taken and stored in timings:

  • evaluation of the Jacobian in jacobian
  • calling the sum fatorization subroutine in whole_sumfact
  • symmetrizing the element matrix in symmetrize_em

◆ timings()

template<class F = Real>
static bool hp3D::Identity< F >::timings ( )
static

Returns true if the class is able to do timings.

The ability to do timings depends on a compiler switch in bilinearform.cc file.

Member Data Documentation

◆ frm_

template<class F = Real>
concepts::ElementFormulaContainer<F> hp3D::Identity< F >::frm_
private

Element formula.

Definition at line 139 of file bilinearForm.hh.

◆ jacobian_

template<class F = Real>
concepts::Array<F> hp3D::Identity< F >::jacobian_
mutableprivate

Intermediate data for element matrix computation.

Definition at line 134 of file bilinearForm.hh.

◆ oldElm_

template<class F = Real>
const Hexahedron* hp3D::Identity< F >::oldElm_
mutableprivate

Intermediate data for element matrix computation.

Definition at line 136 of file bilinearForm.hh.

◆ timeCntr_

template<class F = Real>
uint hp3D::Identity< F >::timeCntr_
staticprivate

Counter for timing table.

Definition at line 144 of file bilinearForm.hh.

◆ timings_

template<class F = Real>
concepts::InOutParameters* hp3D::Identity< F >::timings_
staticprivate

Place to store timing values.

Definition at line 142 of file bilinearForm.hh.


The documentation for this class was generated from the following file:
void addMap(enum mapTypes type, const char *name, const InOutParameters &holder)
Organizes the results in the hashes from InOutParameters in a nice table.
Definition: resultsTable.hh:23
void print(std::ostream &os) const
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich