Skip to content
Snippets Groups Projects
Commit 8a79f686 authored by Naumann, Andreas's avatar Naumann, Andreas
Browse files

added test data creation programs

parent a56db0a5
No related branches found
No related tags found
No related merge requests found
#include "Project.h"
#include "BallProject.h"
#include "io/ArhWriter.h"
#include "io/ArhReader.h"
#include <vector>
using namespace AMDiS;
void write(SolutionInformation& info, std::string filename) {
if(info.sysVec != NULL) {
std::vector< DOFVector< double >* > data(info.sysVec->getNumVectors());
Mesh* mesh = info.sysVec->getDOFVector(0)->getFeSpace()->getMesh();
for(int i=0; i < info.sysVec->getNumVectors(); ++i) {
data[i] = info.sysVec->getDOFVector(i);
}
ArhWriter::write(filename, mesh, data);
}else if(info.dofVec != NULL) {
Mesh* mesh = info.dofVec->getFeSpace()->getMesh();
ArhWriter::write(filename, mesh, info.dofVec);
}else
assert(false);
}
bool compare(SolutionInformation& info, std::string filename) {
if(info.sysVec != NULL) {
assert(false);
}else if(info.dofVec != NULL) {
DOFVector< double > fileVec(info.dofVec->getFeSpace(), "fileVec");
Mesh fileMesh("",2);
fileMesh = *(info.dofVec->getFeSpace()->getMesh());
ArhReader::read(filename, &fileMesh, &fileVec);
}else
assert(false);
return false;
}
#include "SphereProject.h"
#include "AdaptStationary.h"
/// RHS function
class F : public AbstractFunction<double, WorldVector<double> >
{
public:
F(int degree) : AbstractFunction<double, WorldVector<double> >(degree) {}
/// Implementation of AbstractFunction::operator().
double operator()(const WorldVector<double>& x) const
{
return -2.0 * x[0];
}
};
Spheredemo::Spheredemo():
sphere("sphere"),
ballCenter(NULL),
adaptInfo(NULL),
adapt(NULL),
matrixOperator(NULL),
rhsOperator(NULL)
{}
Spheredemo::~Spheredemo() {
if(matrixOperator != NULL)
delete matrixOperator;
if(rhsOperator != NULL)
delete rhsOperator;
if(ballCenter != NULL)
delete ballCenter;
if(adapt != NULL)
delete adapt;
if(adaptInfo != NULL)
delete adaptInfo;
}
void Spheredemo::create(string filename) {
// ===== init parameters =====
Parameters::init(false, filename);
// ===== create projection =====
ballCenter = new WorldVector< double >;
ballCenter->set(0.0);
new BallProject(1, BOUNDARY_PROJECTION, *ballCenter, 1.0);
// ===== create and init the scalar problem =====
sphere.initialize(INIT_ALL);
// === create adapt info ===
adaptInfo = new AdaptInfo("sphere->adapt");
// === create adapt ===
adapt = new AdaptStationary("sphere->adapt",
&sphere,
adaptInfo);
// ===== create matrix operator =====
matrixOperator = new Operator(sphere.getFeSpace());
matrixOperator->addSecondOrderTerm(new Laplace_SOT);
sphere.addMatrixOperator(matrixOperator);
// ===== create rhs operator =====
rhsOperator = new Operator(sphere.getFeSpace());
int degree = sphere.getFeSpace()->getBasisFcts()->getDegree();
rhsOperator->addZeroOrderTerm(new CoordsAtQP_ZOT(new F(degree)));
sphere.addVectorOperator(rhsOperator);
}
int Spheredemo::solve(SolutionInformation& info) {
int retCode = adapt->adapt();
info.dofVec = sphere.getSolution();
return retCode;
}
void createProjectList(ProjectList& list) {
}
#ifndef SPHEREPROJECT_H
#define SPHEREPROJECT_H
#include "Project.h"
#include "ProblemScal.h"
#include "FixVec.h"
class Spheredemo : public Project {
ProblemScal sphere;
WorldVector< double >* ballCenter;
AdaptInfo* adaptInfo;
AdaptStationary* adapt;
Operator* matrixOperator;
Operator* rhsOperator;
public:
Spheredemo();
~Spheredemo();
void create(string filename);
int solve(SolutionInformation& );
};
#endif
#ifndef TORUSPROJECT_H
#define TORUSPROJECT_H
#include "Project.h"
#include "ProblemScal.h"
class Torusdemo : public Project {
ProblemScal torus;
AdaptInfo* adaptInfo;
AdaptStationary* adapt;
Operator* matrixOperator;
Operator* rhsOperator;
public:
Torusdemo();
~Torusdemo();
void create(string filename);
int solve(SolutionInformation&);
};
#endif
#include "VecelliptProject.h"
#include "AdaptStationary.h"
/// Dirichlet boundary function
class G : public AbstractFunction<double, WorldVector<double> >
{
public:
/// Implementation of AbstractFunction::operator().
double operator()(const WorldVector<double>& x) const
{
return exp(-10.0 * (x * x));
}
};
/// RHS function
class F : public AbstractFunction<double, WorldVector<double> >
{
public:
F(int degree) : AbstractFunction<double, WorldVector<double> >(degree) {}
/// Implementation of AbstractFunction::operator().
double operator()(const WorldVector<double>& x) const
{
int dim = x.getSize();
double r2 = (x * x);
double ux = exp(-10.0 * r2);
return -(400.0 * r2 - 20.0 * dim) * ux;
};
};
Vecelliptdemo::Vecelliptdemo():
vecellipt("vecellipt"),
adaptInfo(NULL),
adapt(NULL),
matrixOperator00(NULL),
matrixOperator10(NULL),
matrixOperator11(NULL),
rhsOperator0(NULL)
{}
Vecelliptdemo::~Vecelliptdemo() {
if(matrixOperator00 != NULL)
delete matrixOperator00;
if(matrixOperator10 != NULL)
delete matrixOperator10;
if(matrixOperator11 != NULL)
delete matrixOperator11;
if(rhsOperator0 != NULL)
delete rhsOperator0;
if(adapt != NULL)
delete adapt;
if(adaptInfo != NULL)
delete adaptInfo;
}
void Vecelliptdemo::create(string filename) {
// ===== init parameters =====
Parameters::init(false, filename);
// ===== create and init the scalar problem =====
vecellipt.initialize(INIT_ALL);
// === create adapt info ===
adaptInfo = new AdaptInfo("vecellipt->adapt", vecellipt.getNumComponents());
// === create adapt ===
adapt = new AdaptStationary("vecellipt->adapt", &vecellipt, adaptInfo);
// ===== create matrix operators =====
matrixOperator00 = new Operator(vecellipt.getFeSpace(0), vecellipt.getFeSpace(0));
matrixOperator00->addSecondOrderTerm(new Laplace_SOT);
vecellipt.addMatrixOperator(matrixOperator00, 0, 0);
matrixOperator10 = new Operator(vecellipt.getFeSpace(1), vecellipt.getFeSpace(0));
matrixOperator10->addZeroOrderTerm(new Simple_ZOT);
vecellipt.addMatrixOperator(matrixOperator10, 1, 0);
matrixOperator11 = new Operator(vecellipt.getFeSpace(1), vecellipt.getFeSpace(1));
matrixOperator11->addZeroOrderTerm(new Simple_ZOT(-1.0));
vecellipt.addMatrixOperator(matrixOperator11, 1, 1);
// ===== create rhs operator =====
int degree = vecellipt.getFeSpace(0)->getBasisFcts()->getDegree();
rhsOperator0 = new Operator(vecellipt.getFeSpace(0));
rhsOperator0->addZeroOrderTerm(new CoordsAtQP_ZOT(new F(degree)));
vecellipt.addVectorOperator(rhsOperator0, 0);
// ===== add boundary conditions =====
vecellipt.addDirichletBC(1, 0, 0, new G);
}
int Vecelliptdemo::solve(SolutionInformation& info) {
int retCode = adapt->adapt();
info.sysVec = vecellipt.getSolution();
return retCode;
}
void createProjectList(ProjectList& list) {
}
#ifndef VECELLIPTPROJECT_H
#define VECELLIPTPROJECT_H
#include "Project.h"
#include "ProblemVec.h"
class Vecelliptdemo : public Project {
ProblemVec vecellipt;
AdaptInfo* adaptInfo;
AdaptStationary* adapt;
Operator* matrixOperator00;
Operator* matrixOperator10;
Operator* matrixOperator11;
Operator* rhsOperator0;
public:
Vecelliptdemo();
~Vecelliptdemo();
void create(string filename);
int solve(SolutionInformation& );
};
#endif
#include "Project.h"
#include "@PROJECTINCLUDE@"
int main(int argc, char** argv) {
ProjectList list;
createProjectList(list);
for(ProjectList::iterator it=list.begin() ; it != list.end(); ++it) {
it->create();
SolutionInformation info;
it->solve(info);
write(info, it->getFilename());
}
return 0;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment