Liebe Gitlab-Nutzer, lieber Gitlab-Nutzer, es ist nun möglich sich mittels des ZIH-Logins/LDAP an unserem Dienst anzumelden. Ein Anmelden über dieses erzeugt ein neues Konto. Das alte Konto ist über den Reiter "Standard" erreichbar. Die Administratoren

Dear Gitlab user, it is now possible to log in to our service using the ZIH login/LDAP. Logging in via this will create a new account. The old account can be accessed via the "Standard" tab. The administrators

Commit 3119b7a7 authored by Thomas Witkowski's avatar Thomas Witkowski
Browse files

Updated Reinit library to some changes in AMDiS.

parent d2b344e6
......@@ -9,21 +9,17 @@
using namespace AMDiS;
class BoundaryElementDist {
public:
MEMORY_MANAGED(BoundaryElementDist);
public:
/**
* Constructor.
*/
BoundaryElementDist(ElementLevelSet *elLS_, int dim_)
: dim(dim_),
elLS(elLS_)
elLS(elLS_)
{
FUNCNAME("BoundaryElementDist::BoundaryElementDist()");
TEST_EXIT(dim == 2 || dim == 3)
("function only works for dimension 2 !\n");
};
}
/**
* Virtual destructor.
......
......@@ -12,16 +12,12 @@ using namespace AMDiS;
class BoundaryElementEdgeDist : public BoundaryElementDist
{
public:
MEMORY_MANAGED(BoundaryElementEdgeDist);
public:
/**
* Constructor.
*/
BoundaryElementEdgeDist(ElementLevelSet *elLS_,
int dim_)
: BoundaryElementDist(elLS_, dim_)
{};
{}
/**
* Calculates distance from the interface for all vertices of a boundary
......
......@@ -3,24 +3,18 @@
#include "ElInfo.h"
#include "FixVec.h"
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
using namespace AMDiS;
class BoundaryElementLevelSetDist : public BoundaryElementDist
{
public:
MEMORY_MANAGED(BoundaryElementLevelSetDist);
public:
/**
* Constructor.
*/
BoundaryElementLevelSetDist(ElementLevelSet *elLS_, int dim_)
: BoundaryElementDist(elLS_, dim_)
{};
{}
/**
* Calculates distance from the interface for all vertices of a boundary
......
......@@ -3,24 +3,18 @@
#include "ElInfo.h"
#include "FixVec.h"
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
using namespace AMDiS;
class BoundaryElementNormalDist : public BoundaryElementDist
{
public:
MEMORY_MANAGED(BoundaryElementNormalDist);
public:
/**
* Constructor.
*/
BoundaryElementNormalDist(ElementLevelSet *elLS_, int dim_)
: BoundaryElementDist(elLS_, dim_)
{};
{}
/**
* Calculates distance from the interface for all vertices of a boundary
......
......@@ -3,9 +3,7 @@
#include "ElInfo.h"
#include "FixVec.h"
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
#include "VelocityExt.h"
......@@ -13,24 +11,16 @@ using namespace AMDiS;
class BoundaryElementTopDist : public BoundaryElementDist
{
public:
MEMORY_MANAGED(BoundaryElementTopDist);
/**
* Constructor.
*/
public:
BoundaryElementTopDist(ElementLevelSet *elLS_,
int dim_,
VelocityExt *velExt_ = NULL)
: BoundaryElementDist(elLS_, dim_),
velExt(velExt_)
{};
velExt(velExt_)
{}
/**
* Destructor.
*/
~BoundaryElementTopDist()
{};
{}
/**
* Calculates distance from the interface for all vertices of a boundary
......
#ifndef AMDIS_ELEMENTLEVELSET_H
#define AMDIS_ELEMENTLEVELSET_H
#include "AMDiS_fwd.h"
#include "ElementFunction.h"
#include "FixVec.h"
// #include "MemoryManager.h"
#include "Parameters.h"
namespace AMDiS {
class Element;
class ElInfo;
class Mesh;
}
using namespace AMDiS;
using namespace std;
......@@ -58,26 +52,21 @@ using namespace std;
class ElementLevelSet
{
public:
MEMORY_MANAGED(ElementLevelSet);
/**
* Constructor.
*/
ElementLevelSet(const char *name_,
ElementFunction<double> *lSFct_,
Mesh *mesh_)
: name(name_),
elInfo(NULL),
lastEl(NULL),
level_set_domain(LEVEL_SET_UNDEFINED),
numIntersecPoints(0),
elStatus(LEVEL_SET_UNDEFINED),
numElVertexInterior(0),
numElVertexBoundary(0),
numElVertexExterior(0),
LS_VAL_TOL(1.e-8),
LS_VAL_MIN(1.e-8),
SP_BARY_TOL(1.e-7)
elInfo(NULL),
lastEl(NULL),
level_set_domain(LEVEL_SET_UNDEFINED),
numIntersecPoints(0),
elStatus(LEVEL_SET_UNDEFINED),
numElVertexInterior(0),
numElVertexBoundary(0),
numElVertexExterior(0),
LS_VAL_TOL(1.e-8),
LS_VAL_MIN(1.e-8),
SP_BARY_TOL(1.e-7)
{
FUNCNAME("ElementLevelSet::ElementLevelSet()");
......@@ -108,11 +97,8 @@ class ElementLevelSet
TEST_EXIT(LS_VAL_MIN > 0)("illegal LS_VAL_MIN\n");
TEST_EXIT(SP_BARY_TOL > 0)("illegal SP_BARY_TOL\n");
}
};
}
/**
* Destructor.
*/
~ElementLevelSet()
{
if (elVertexStatusVec)
......@@ -121,7 +107,7 @@ class ElementLevelSet
delete [] elVertexLevelSetVec;
if (elIntersecPoints)
DELETE elIntersecPoints;
};
}
/**
* Calculates LevelSet-status of element and its intersection points
......
......@@ -2,28 +2,18 @@
#define ELEMENTUPDATE_H
#include "FixVec.h"
// #include "MemoryManager.h"
#include "VelocityExt.h"
using namespace AMDiS;
class ElementUpdate
{
public:
MEMORY_MANAGED(ElementUpdate);
/**
* Constructor.
*/
public:
ElementUpdate(VelocityExt *velExt_)
: velExt(velExt_)
{};
{}
/**
* Virtual destructor.
*/
virtual ~ElementUpdate() {};
virtual ~ElementUpdate() {}
/**
* Pure virtual function.
......
......@@ -2,8 +2,6 @@
#define ELEMENTUPDATE_2D_H
#include "FixVec.h"
// #include "MemoryManager.h"
#include "ElementUpdate.h"
#include "VelocityExt.h"
......@@ -11,15 +9,10 @@ using namespace AMDiS;
class ElementUpdate_2d : public ElementUpdate
{
public:
MEMORY_MANAGED(ElementUpdate_2d);
/**
* Constructor.
*/
public:
ElementUpdate_2d(VelocityExt *velExt_ = NULL)
: ElementUpdate(velExt_)
{};
{}
/**
* Realization of ElementUpdate::calcElementUpdate.
......
......@@ -2,8 +2,6 @@
#define ELEMENTUPDATE_3D_H
#include "FixVec.h"
// #include "MemoryManager.h"
#include "ElementUpdate.h"
#include "ElementUpdate_2d.h"
#include "VelocityExt.h"
......@@ -12,25 +10,17 @@ using namespace AMDiS;
class ElementUpdate_3d : public ElementUpdate
{
public:
MEMORY_MANAGED(ElementUpdate_3d);
/**
* Constructor.
*/
public:
ElementUpdate_3d(VelocityExt *velExt_ = NULL)
: ElementUpdate(velExt_)
{
elUpdate2d = NEW ElementUpdate_2d(velExt_);
};
elUpdate2d = new ElementUpdate_2d(velExt_);
}
/**
* Destructor.
*/
~ElementUpdate_3d()
{
DELETE elUpdate2d;
};
delete elUpdate2d;
}
/**
* Realization of ElementUpdate::calcElementUpdate.
......
......@@ -8,11 +8,8 @@
#include "FileWriter.h"
#include "FixVec.h"
#include "Flag.h"
// #include "MemoryManager.h"
#include "Parameters.h"
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
#include "BoundaryElementLevelSetDist.h"
#include "BoundaryElementTopDist.h"
......@@ -41,12 +38,7 @@ using namespace AMDiS;
//////////////////////////////////////////////////////////////////////////////
class HL_SignedDist
{
public:
MEMORY_MANAGED(HL_SignedDist);
/**
* Constructor.
*/
public:
HL_SignedDist(const char *name_,
int dim_,
bool doVelocityExt = false,
......@@ -77,11 +69,11 @@ class HL_SignedDist
// ===== Create functionality for velocity extension. =====
if (doVelocityExt) {
if (velExtType.isSet(VEL_EXT))
velExt = NEW VelocityExt(dim);
velExt = new VelocityExt(dim);
else
velExt = NEW VelocityExtFromVelocityField(dim);
velExt = new VelocityExtFromVelocityField(dim);
}
};
}
/**
* Virtual destructor.
......@@ -89,19 +81,19 @@ class HL_SignedDist
virtual ~HL_SignedDist()
{
if (elUpdate)
DELETE elUpdate;
delete elUpdate;
if (bndElDist)
DELETE bndElDist;
delete bndElDist;
if (elLS)
DELETE elLS;
delete elLS;
if (phi)
DELETE phi;
delete phi;
if (bound_DOF)
DELETE bound_DOF;
delete bound_DOF;
DELETE velExt;
delete velExt;
};
/**
......@@ -251,13 +243,13 @@ class HL_SignedDist
*/
void printLevelSetFct()
{
FileWriter *fileWriter = NEW FileWriter(
FileWriter *fileWriter = new FileWriter(
"SignedDist->level set fct output",
feSpace->getMesh(),
const_cast<DOFVector<double> *>(lS_DOF));
fileWriter->writeFiles(adaptInfo, false);
DELETE fileWriter;
delete fileWriter;
};
/**
......@@ -265,13 +257,13 @@ class HL_SignedDist
*/
void printSignedDistFct()
{
FileWriter *fileWriter = NEW FileWriter(
FileWriter *fileWriter = new FileWriter(
"SignedDist->result output",
feSpace->getMesh(),
sD_DOF);
fileWriter->writeFiles(adaptInfo, false);
DELETE fileWriter;
delete fileWriter;
};
protected:
......@@ -306,13 +298,12 @@ class HL_SignedDist
*/
void printBoundInitFct()
{
FileWriter *fileWriter = NEW FileWriter(
"SignedDist->boundary initialization output",
feSpace->getMesh(),
sD_DOF);
FileWriter *fileWriter = new FileWriter("SignedDist->boundary initialization output",
feSpace->getMesh(),
sD_DOF);
fileWriter->writeFiles(adaptInfo, false);
DELETE fileWriter;
delete fileWriter;
};
public:
......
#ifndef HL_SIGNEDDISTBORNEMANN
#define HL_SIGNEDDISTBORNEMANN
#include <queue>
#include <time.h>
#include "ElInfo.h"
#include "FixVec.h"
// #include "MemoryManager.h"
#include "Traverse.h"
#include <queue>
#include <time.h>
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
/* #include "BoundaryElementTopDist.h" */
/* #include "BoundaryElementEdgeDist.h" */
/* #include "BoundaryElementNormalDist.h" */
#include "ElementUpdate.h"
#include "ElementUpdate_2d.h"
#include "ElementUpdate_3d.h"
......@@ -32,25 +26,20 @@ typedef struct
class HL_SignedDistBornemann : public HL_SignedDist
{
public:
MEMORY_MANAGED(HL_SignedDistBornemann);
/**
* Constructor.
*/
public:
HL_SignedDistBornemann(const char *name_,int dim_)
: HL_SignedDist(name_, dim_),
smiAdapter(NULL)
smiAdapter(NULL)
{
FUNCNAME("HL_SignedDistBornemann::HL_SignedDistBornemann");
// ===== Read parameters from init file. =====
GET_PARAMETER(0,name + "->tolerance", "%f", &tol);
GET_PARAMETER(0,name + "->count_how_often_saved_in_list", "%d", &count_in_list);
GET_PARAMETER(0,name + "->save_in_list->the ..th", "%d", &print_in_list);
GET_PARAMETER(0, name + "->tolerance", "%f", &tol);
GET_PARAMETER(0, name + "->count_how_often_saved_in_list", "%d", &count_in_list);
GET_PARAMETER(0, name + "->save_in_list->the ..th", "%d", &print_in_list);
TEST_EXIT(tol > 0)("illegal tolerance !\n");
};
}
protected:
/**
......
#ifndef HL_SIGNEDDISTLEVELS
#define HL_SIGNEDDISTLEVELS
#include "ElInfo.h"
#include "FixVec.h"
// #include "MemoryManager.h"
#include "Traverse.h"
#include <queue>
#include <vector>
#include <time.h>
#include "ElInfo.h"
#include "FixVec.h"
#include "Traverse.h"
#include "ElementLevelSet.h"
#include "BoundaryElementDist.h"
#include "ElementUpdate.h"
#include "ElementUpdate_2d.h"
#include "ElementUpdate_3d.h"
#include "HL_SignedDist.h"
#include "VelocityExt.h"
#include "SMIAdapter.h"
#include "smi.h"
......@@ -31,31 +27,26 @@ typedef struct
class HL_SignedDistLevels : public HL_SignedDist
{
public:
MEMORY_MANAGED(HL_SignedDistLevels);
/**
* Constructor.
*/
public:
HL_SignedDistLevels(const char *name_,
int dim_,
bool doVelocityExt = false,
Flag velExtType_ = VEL_EXT)
: HL_SignedDist(name_, dim_, doVelocityExt, velExtType_),
smiAdapter(NULL)
smiAdapter(NULL)
{
FUNCNAME("HL_SignedDistLevels::HL_SignedDistLevels");
FUNCNAME("HL_SignedDistLevels::HL_SignedDistLevels()");
// ===== Read parameters from init file. =====
GET_PARAMETER(0,name + "->tolerance", "%f", &tol);
GET_PARAMETER(0,name + "->count_how_often_saved_in_list", "%d", &count_in_list);
GET_PARAMETER(0,name + "->save_in_list->the ..th", "%d", &print_in_list);
GET_PARAMETER(0,name + "->save_in_list->after the ..th traversing of the list", "%d", &print_in_list_2);
GET_PARAMETER(0,name + "->print_level", "%d", &chosen_level);
GET_PARAMETER(0,name + "->print_level_yes_no", "%d", &print_level);
GET_PARAMETER(0, name + "->tolerance", "%f", &tol);
GET_PARAMETER(0, name + "->count_how_often_saved_in_list", "%d", &count_in_list);
GET_PARAMETER(0, name + "->save_in_list->the ..th", "%d", &print_in_list);
GET_PARAMETER(0, name + "->save_in_list->after the ..th traversing of the list", "%d", &print_in_list_2);
GET_PARAMETER(0, name + "->print_level", "%d", &chosen_level);
GET_PARAMETER(0, name + "->print_level_yes_no", "%d", &print_level);
TEST_EXIT(tol > 0)("illegal tolerance !\n");
};
}
protected:
/**
......
#include "HL_SignedDistTraverse.h"
#include "VelocityExtFromVelocityField.h"
void
HL_SignedDistTraverse::initializeBoundary()
void HL_SignedDistTraverse::initializeBoundary()
{
FUNCNAME("HL_SignedDistTraverse::initializeBoundary()");
......@@ -16,7 +14,7 @@ HL_SignedDistTraverse::initializeBoundary()
int elStatus;
const int nBasFcts = feSpace->getBasisFcts()->getNumber();
DegreeOfFreedom *locInd = GET_MEMORY(DegreeOfFreedom, nBasFcts);
DegreeOfFreedom *locInd = new DegreeOfFreedom[nBasFcts];
ElInfo *elInfo;
if (velExt && velExtType.isSet(VEL_EXT_FROM_VEL_FIELD)) {
......@@ -26,8 +24,7 @@ HL_SignedDistTraverse::initializeBoundary()
Mesh::FILL_BOUND |
Mesh::FILL_COORDS |
Mesh::FILL_GRD_LAMBDA);
}
else {
} else {
elInfo = stack.traverseFirst(feSpace->getMesh(),
-1,
Mesh::CALL_LEAF_EL |
......@@ -94,9 +91,7 @@ HL_SignedDistTraverse::initializeBoundary()
elInfo = stack.traverseNext(elInfo);
} // end of: mesh traverse
FREE_MEMORY(locInd, DegreeOfFreedom, nBasFcts);
return;
delete [] locInd;
}
void
......@@ -165,7 +160,7 @@ HL_SignedDistTraverse::HL_elementUpdate(ElInfo *elInfo)
// ===== Get global indices of vertices of element. =====
const int nBasFcts = feSpace->getBasisFcts()->getNumber();
DegreeOfFreedom *locInd = GET_MEMORY(DegreeOfFreedom, nBasFcts);
DegreeOfFreedom *locInd = new DegreeOfFreedom[nBasFcts];
feSpace->getBasisFcts()->getLocalIndices(
const_cast<Element *>(elInfo->getElement()),
......@@ -203,13 +198,12 @@ HL_SignedDistTraverse::HL_elementUpdate(ElInfo *elInfo)
}
}
FREE_MEMORY(locInd, DegreeOfFreedom, nBasFcts);
delete [] locInd;
}
double
HL_SignedDistTraverse::calcElementUpdate(ElInfo *elInfo,
int nXh,
const DegreeOfFreedom *locInd)
double HL_SignedDistTraverse::calcElementUpdate(ElInfo *elInfo,
int nXh,
const DegreeOfFreedom *locInd)
{
double update;
FixVec<WorldVector<double> *, VERTEX> elVert(dim, NO_INIT);
......@@ -257,8 +251,7 @@ HL_SignedDistTraverse::calcElementUpdate(ElInfo *elInfo,
return update;
}
bool
HL_SignedDistTraverse::checkTol()
bool HL_SignedDistTraverse::checkTol()