Skip to content
Snippets Groups Projects
Commit ede0268d authored by Oliver Sander's avatar Oliver Sander Committed by sander
Browse files

Get the initial deformation from a python lambda

[[Imported from SVN: r9797]]
parent e696e4a1
No related branches found
No related tags found
No related merge requests found
......@@ -88,4 +88,7 @@ kappa = 1
# x is the vertex coordinate
dirichletVerticesPredicate = "x[0] < 0.001 or x[0] > 0.0999"
# Initial deformation
initialDeformation = "[x[0], x[1], 0]"
......@@ -137,3 +137,6 @@ neumannValues = -1e3 0 0
structuredGrid = false
path = /home/sander/data/shells/wriggers_L_shape/
gridFile = wriggers-L-shape_99.msh
# Initial deformation
initialDeformation = "[x[0], x[1], 0.002*cos(1e4*x[0])]"
......@@ -95,3 +95,5 @@ neumannVerticesPredicate = "x[1] < -0.239"
### Neumann values, if needed
neumannValues = -1e3 0 0
# Initial deformation
initialDeformation = "[x[0], x[1], 0]"
......@@ -51,21 +51,6 @@ const int blocksize = TargetSpace::TangentVector::dimension;
using namespace Dune;
class Identity
: public Dune::VirtualFunction<FieldVector<double,dim>, FieldVector<double,3>>
{
public:
void evaluate(const FieldVector<double,dim>& x, FieldVector<double,3>& y) const
{
y = 0;
for (int i=0; i<dim; i++)
y[i] = x[i];
//y[2] = 0.002*std::cos(1e4*x[0]);
}
};
// Dirichlet boundary data for the shear/wrinkling example
class WrinklingDirichletValues
: public Dune::VirtualFunction<FieldVector<double,dim>, FieldVector<double,3> >
......@@ -312,9 +297,11 @@ int main (int argc, char *argv[]) try
SolutionType x(feBasis.size());
Identity identity;
lambda = std::string("lambda x: (") + parameterSet.get<std::string>("initialDeformation") + std::string(")");
PythonFunction<FieldVector<double,dim>, FieldVector<double,3> > pythonInitialDeformation(Python::evaluate(lambda));
std::vector<FieldVector<double,3> > v;
Functions::interpolate(feBasis, v, identity);
Functions::interpolate(feBasis, v, pythonInitialDeformation);
for (size_t i=0; i<x.size(); i++)
x[i].r = v[i];
......
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