Department of Scientific Computing   
Institute for Numerical Simulation   
University of Bonn   
Documentation
Download
Programming References
Bug Reports / Suggestions
FAQ
Authors
Main Page   Class Hierarchy   Alphabetical List   Compound List   File List   Compound Members  

AdaptiveData Struct Template Reference

AdaptiveData<DIM> is just a collection of ordinary vectors of doubles which contain the actual numerical data associated to an AdaptiveGrid<DIM>. More...

#include <AdaptiveData.hpp>

List of all members.

Public Methods

 AdaptiveData ()
 default constructor.

 AdaptiveData (AdaptiveGrid< DIM > *G)
 constructor with automatic attachment to an adaptive grid.

void Attach (AdaptiveGrid< DIM > *G)
 attach to a certain grid.

void DeAttach ()
 de-attach.

size_t Size ()
 return number of entries.

size_t AllocatedEntries ()
 return number of allocated entries.

void SetRefine ()
 The data must be attached to a grid, therefore, refinement of this grid affects the data.
If SetRefine() has been called, grid refinement automatically leads to a 'refinement' of the data: coefficients for indices which are in the old and the new grid are retained, coefficients for completely new indices are set to zero.


void FromUniform (UniformData< DIM > *A, int *J)
 initialize AdaptiveGrid connected to (*this) as uniform grid up to level J[.] using AdaptiveGrid<DIM>::SetFull(), and copy data to (*this).

void ToUniform (UniformData< DIM > *A, int *J)
 write (*this) to a uniform data structure.

double* GetP (int *l, int *t)
 Return pointer to coefficient (l,t), return NULL if illegal index.

int Set (int *l, int *t, const double a)
 coefficient(l,t)==a, return 0 for successfull op / !=0 for error code.

double Get (int *l, int *t)
 Return coefficient (l,t).

void Set (int *l, const double a)
 Set all coefficients of level *l to a.

void Set (int dir, int k, int s, const double a)
 Set all coefficients (l,t) whith l[dir]==k and t[dir]==s to a.

void Get (int *l, int *t, int dir, AdaptiveB *B, int leftrightall=0)
 write the coefficients of a 1-dimensional line through (l,t) along the dir-th coordinate direction in *B
if leftrightall=0 then all coefficients are written to *B
if ==1 then just those from the left boundary
if ==2 then just those from the right boundary
if ==3 then just those from both boundaries.


void Set (int *l, int *t, int dir, AdaptiveB *B, int leftrightall=0, bool add=false)
 read the coefficients of a 1-dimensional line through (l,t) along the dir-th coordinate direction from *B if add==true, add the coefficients from *B.

void ReadSlice (AdaptiveData< DIM+1 > *X)
 *this must be attached to an AdaptiveGrid which is a slice of the (DIM+1)-dimensional AdaptiveGrid *X is attached to. Coefficients from *X are copied to *this.

void WriteSlice (AdaptiveData< DIM+1 > *X, bool add=false)
 *this must be attached to an AdaptiveGrid which is a slice of the (DIM+1)-dimensional AdaptiveGrid *X is attached to. If add==true coefficients from *this are added to *X; otherwise they are copied to *X.

void SetFunction (Function *F, bool boundaryonly=false, bool notransform=false)
 Compute the interpolant of (*F) with respect to the given grid. Supported for Interpolets only. On periodic faces (known because of the link to AdaptiveGrid!) the interpolant has periodic BC; on all other faces the boundary conditions are set BC[][]=-1, i.e. the underlying basis functions have no BC there.
If (boundaryonly==true) then a quite smooth extension of JUST the boundary values is computed. If notransform==false the nodal values are wavelet transformed, otherwise we keep the nodal represenation. In the later case, boundaryonly must be false.


void SetBoundaryValueFunction (AdaptiveData< DIM-1 > *F[DIM][2], int BC[DIM][2], bool ContinuousDirichletValues=false)
 Compute a function which approximately takes the given Dirichlet- or Neumann-boundary conditions. Supported for Interpolets only. The result has either no (BC[][]=-1) or periodic boundary conditions on the respective faces.
In case of inconsistencies of Dirichlet-values on common (DIM-2)-dimensional faces of neighbouring Dirichlet-faces some avereging is performed by default.
A similar problem occurs, if a Dirichlet- and a Neumann face are neighbours, but the tangential derivatives of the Dirichlet values do not match the Neumann values on the common (DIM-2)-dimensional face. In this case we assume that Dirichlet values have a higher priority, and the Neumann values are altered appropriately in the closest possible vicinity of the common face.

Attention: the current implementation can handle corners only, where at most three Neumann-faces touch. (The number of also touching Dirichlet faces doesn't matter.)
Parameters:
F   are pointers to either Dirichlet or Neumann boundary values, F[i][k] must be attached to G->FaceGrid[i][k]. In contrast to the usual convention, Neumann conditions do not mean the normal derivative, but just the derivative on the resptive faces, i.e. there is no change of sign e.g. at the face $x=0$ and the face $x=1$. This convention is more convenient for our purposes and often saves the rather artificial discrimination between left and right faces.
BC   determines which boundary conditions (BC[][]==0 means Dirichlet, BC[][]==1 Neumann, BC[][]=-1 means don't care about this face)
ContinuousDirichletValues   If two (DIM-1)-dimensional Dirichlet faces touch, the Dirichlet boundary values prescribed there may not collapse for both faces. The Dirichlet values are not continuous then. If ContinuousDirichletValues == false, extra effort is made to enforce collapsing boundary values by means of averaging nodal values of the Dirichlet values.
This involves wavelet transforms (and inverse transforms) of the boundary values. The Dirichlet values may be given either as wavelet coefficients or nodal values. SetBoundaryValueFunction automatically generates the representation required internally. E.g., if ContinuousDirichletValues=false it is advantageous to give the boundary values in nodal representation, as some extra wavelet transforms may be skipped then, which saves computing time.


double Max ()
 return maximum value.

double Min ()
 return minimum value.

double MaxAbs ()
 return maximum absolute value.

double Sum ()
 return sum of all entries.

double SumAbs ()
 return sum of absolute entries.

double SumSqr ()
 return sum of all squared entries.

double LpNorm (AdaptiveData< DIM > *Tmp, double p)
 calculate (approximate) Lp-Norms via transform in physical space, pointwise evaluation of |.|^p transform back in wavelet space and exact integration of this interpolant
requires auxiliary data structure, *Tmp==*this is valid, but then *this is altered
use p==HUGE_VAL to get maximum norm.


double LmaxNorm (AdaptiveData< DIM > *Tmp)
 maximum norm.

double L2Norm (AdaptiveData< DIM > *Tmp)
 L^2 norm.

double L1Norm (AdaptiveData< DIM > *Tmp)
 L^1 norm.

double Integrate ()
 compute exact integral of *this.

double InnerProd (AdaptiveData< DIM > *X)
 compute inner product of the entries of *this and *X.

void Set (const double a)
 all entries are set to a.

void Copy (AdaptiveData< DIM > *X)
 copy everything from *X.

void Abs (AdaptiveData< DIM > *X)
 compute componentwise the absolute values.

void Pow (AdaptiveData< DIM > *X, double p)
 compute componentwise the given power.

void Sqr (AdaptiveData< DIM > *X)
 compute componentwise the square ( faster than pow(X,2) ).

void PPlus (AdaptiveData< DIM > *X)
 (*this) += (*X).

void PPlus (const double a, AdaptiveData< DIM > *X)
 (*this) += a*(*X).

void Add (AdaptiveData< DIM > *X, AdaptiveData< DIM > *Y)
 (*this) = (*X)+(*Y).

void Add (const double a, AdaptiveData< DIM > *X, const double b, AdaptiveData< DIM > *Y)
 (*this) = a*(*X) + b*(*Y).

void Add (const double a, AdaptiveData< DIM > *X, const double b, AdaptiveData< DIM > *Y, const double c, AdaptiveData< DIM > *Z)
 (*this) = a*(*X) + b*(*Y) + c*(Z).

void Add (const double a, AdaptiveData< DIM > *X, const double b, AdaptiveData< DIM > *Y, const double c, AdaptiveData< DIM > *Z, const double d, AdaptiveData< DIM > *Q)
 (*this) = a*(*X) + b*(*Y) + c*(Z) +d*(Q).

void Add (const double a, AdaptiveData< DIM > *X, const double b, AdaptiveData< DIM > *Y, const double c, AdaptiveData< DIM > *Z, const double d, AdaptiveData< DIM > *Q, const double f, AdaptiveData< DIM > *R)
 (*this) = a*(*X) + b*(*Y) + c*(Z) +d*(Q) +e*(R).

void MMinus (AdaptiveData< DIM > *X)
 (*this) -= (*X).

void MMinus (const double a, AdaptiveData< DIM > *X)
 (*this) -= a*(*X).

void Sub (AdaptiveData< DIM > *X, AdaptiveData< DIM > *Y)
 (*this) = (*X) - (*Y).

void Mul (AdaptiveData< DIM > *X, AdaptiveData< DIM > *Y)
 (*this) = (*X) * (*Y) (componentwise).

void Mul (AdaptiveData< DIM > *X, AdaptiveData< DIM > *Y, double f)
 (*this) = (*X) * (*Y) *f (componentwise).

void Mul (AdaptiveData< DIM > *X, double f)
 (*this) = f (*X).

void Mul (double f)
 (*this) *= f.

void ApplyOp (int *BCT, AdaptiveData< DIM > *X, int dir, unsigned int op)
 apply operation to *X with respect to one coordinate direction
Parameters:
X   argument for the operation; may be *this
BCT   boundary conditions of test functions along the <dir>th coordinate direction and also the boundary conditions of *this after the operation w.r.t. <dir>. The resaon is that e.g. FIRTSDERIVATIVE and also finite difference operations are considered as projections of something onto a test space
dir   coordinate direction with respect to which the operation is applied
op   identifier of the univariate operation valid IDs are:
NOOPERATION do nothing
TRANSFORM wavelet transform, supported for Interpolets only!
INVERSETRANSFORM inverse transform, supported for Interpolets only!
INTERPOLET2LIFTING change of basis from Interpolets to Lifting-Interpolets
LIFTING2INTERPOLET change of basis from
FIRSTDERIVATIVE compute (Petrov-)Galerkin projection of first derivative,
test functions are the dual wavelets, a finite difference approximation is used for low order Interpolets
STIFFNESS as above, but second derivative
PROJECTION as above, just projection
FD12 finite fifferences 1. derivative 2. order accuracy
FD14 finite fifferences 1. derivative 4. order accuracy
FD16 finite fifferences 1. derivative 6. order accuracy
FD22 finite fifferences 2. derivative 2. order accuracy
FD24 finite fifferences 2. derivative 4. order accuracy
FD40 finite fifferences 4. derivative 0. order accuracy
FD60 finite fifferences 6. derivative 0. order accuracy
FD12II as above, but consistent treatment at interval bounds
FD22II as above
FD14II as above
FD24II as above
DIV stabilized div operator</td>
GRAD stabilized div operator</td>
DIVGRAD stabilized div-grad operator</td>
DIVFD4 stabilized div operator with 4th order f.d. stencil
DIVFD6 stabilized div operator with 6th order f.d. stencil
GRADFD4 as above
GRADFD6 as above
DIVGRADFD6 as above
DIVGRADFD4NOSTAB no stabilization
WENO3 3rd order WENO finite difference scheme
WENO5 5th order WENO finite difference scheme
.


void ApplyOp (AdaptiveData< DIM > *X, int dir, unsigned int op)
 as above, but keep boundary conditions of X, i.e. short cut for ApplyOp(X->Ext.BC[dir] , X, dir , op) ;.

void WriteUDF (const char *name, UniformData< DIM > *M, AdaptiveData< DIM > *Tmp=NULL, bool CastToFloat=false)
 Write a file which can be read from VTK (vtkRectilinearGridReader) or in MATLAB using ReadUDF; do help ReadUDF for more MATLAB-information
Parameters:
Tmp   if !=NULL nodal values are written: (*this) is projected to general boundary conditions first, then copied to *M and inverse transformed Tmp may be the same as this
CastToFloat   if true just write float values to save disk space.


void WriteUDF (const char *name, int *L, AdaptiveData< DIM > *Tmp=NULL, bool CastToFloat=false)
 as above, but internal allocation of *M with size 1<<L[i] along each coordinate direction.

void WriteUDF (const char *name, int L, AdaptiveData< DIM > *Tmp=NULL, bool CastToFloat=false)
 as above, but level L for all coordinate directions.

void WriteSparse (const char *name, AdaptiveData< DIM > **M, int N, bool CastToFloat=false, AdaptivityCriterion *R=NULL)
 WriteSparse writes *M[0],..,*M[N-1] in one file. M[0...N-1] must be attached to the same grid This multiple storage saves a lot of disk space as the index information needs to be saved only once. Additionally one can use CastToFloat=true, then the coefficients are saved as floats. If R!=NULL just those coefficients are written which fulfill the refinement criterion.

void WriteSparse (const char *name, AdaptivityCriterion *R=NULL)
 just write *this.

void ReadSparse (const char *name, AdaptiveData< DIM > **M, int N, int *which)
 ReadSparse reads from a file the components which[0],...,which[N-1] and stores them in M[0],...,M[N-1]. M[0],...,M[N-1] must be linked to the same grid and this grid is refined accordingly. You MUST call AdaptiveGrid::ReadSparse() first, to initialize the AdaptiveGrid!!!

void ReadSparse (const char *name)
 read just one Data from file.

void SetBoundaryConditions (AdaptiveData< DIM > *X)
 Copy boundary conditions from *X.

void SetBoundaryConditions (int bc[DIM][2])
 Set boundary conditions.

bool SameBoundaryConditions (AdaptiveData< DIM > *X)
 Return true if *this and *x have the same boundary conditions.

void CopyExtensions (AdaptiveData< DIM > *X)
 copy just the Extensions.

template<int DIM> struct AdaptiveData


The documentation for this struct was generated from the following files:
Generated at Mon Aug 19 10:02:32 2002 for AWFD by doxygen1.2.8.1 written by Dimitri van Heesch, © 1997-2001