← Back to team overview

dolfin team mailing list archive

[Question #157032]: slow interpolation of expression

 

New question #157032 on DOLFIN:
https://answers.launchpad.net/dolfin/+question/157032

Hi @all

I use a semi-lagrange method for some transport model. With this method, I have to evaluate a function at some arbitrary point. It seems to me, that this interpolation / evaluation is really slow, so I wrote a small program for timing.

In this program I measure the time of the  following operations:
     - assign a constant to a function
     - solve the variational problem inner(u,v)*dx = inner(f+g, v)*dx
     - assigh the expression u = f + 2.0
     - assign the expression u = f + g

The variational problem is the ufl code (I called the file Poisson.ufl, that's why I include Poisson.h later):

element = FiniteElement("Lagrange", triangle, 1)

u = TrialFunction(element)
v = TestFunction(element)
f = Coefficient(element)
g = Coefficient(element)

a = inner(u, v)*dx
L = (f+g)*v*dx 

and the c++ code for testing is:

#include <dolfin.h>
#include <vector>
#include <iostream>
#include <time.h>

#include "ufl/Poisson.h"

using namespace dolfin;

typedef Array< double > A;

struct Add : Expression {
	Function& f;
	Add(Function& f):
		Expression(), f(f)
	{}

	void eval(A& values, const A& x) const
	{
		f.eval(values, x);
		for(unsigned int i=0; i < values.size(); ++i)
			values[i] += 2.0;
	}
};

struct EvalSome : Expression {
	Function& f;
	mutable A movedX;
	EvalSome(Function& f):
		f(f)
	{}

	void eval(A& values, const A& x) const 
	{
		//TODO: boundary
		//for(unsigned int i(0); i < x.size(); ++i)
		//	movedX[i] = x[i] - 1.0;
		//f.eval(values, movedX);
		f.eval(values, x);
	}
};

struct Set2 : Expression {
	void eval(A& v, const A&) const 
	{
		v[0] = 1.0;
	}
};

int main(int argc, char* argv[]) {
	Mesh mesh(Rectangle(0.0, 0.0, 1.0, 2.0, 20, 50));
	Poisson::FunctionSpace space(mesh);
	Poisson::BilinearForm a(space, space);
	Poisson::LinearForm l(space);

	Function u(space);
	Function u2(space);

	Function f(space);
	Function g(space);
	l.f = f;
	l.g = g;

	Constant zero(0.0);
	Set2 set2;
	f = zero;
	g = set2;

	VariationalProblem p(a, l);

	timespec beg;
	timespec end;
	clock_gettime(CLOCK_REALTIME, &beg);
	p.solve(u);
	clock_gettime(CLOCK_REALTIME, &end);
	std::cout << "needed " << end.tv_nsec - beg.tv_nsec << " nanoseconds for solve\n";

	Add add(u);
	clock_gettime(CLOCK_REALTIME, &beg);
	u2 = add;
	clock_gettime(CLOCK_REALTIME, &end);
	std::cout << "needed " << end.tv_nsec - beg.tv_nsec << " nanoseconds for setting a dependend expression\n";

	clock_gettime(CLOCK_REALTIME, &beg);
	u2 = set2;
	clock_gettime(CLOCK_REALTIME, &end);
	std::cout << "needed " << end.tv_nsec - beg.tv_nsec << " nanoseconds for setting a constant expression\n";

	clock_gettime(CLOCK_REALTIME, &beg);
	u2 = f;
	clock_gettime(CLOCK_REALTIME, &end);
	std::cout << "needed " << end.tv_nsec - beg.tv_nsec << " nanoseconds for setting an constant\n";

	GenericVector& uVec = u.vector();
	clock_gettime(CLOCK_REALTIME, &beg);
	uVec -= u2.vector();
	clock_gettime(CLOCK_REALTIME, &end);
	std::cout << "needed " << end.tv_nsec - beg.tv_nsec << " nanoseconds for vec sub\n";
	return 0;
}

When I run this program, it showed me, that the command u2 = add, (i.e. the assignment of the Add-struct to the function ) is much slower than the even assembly and solve the variational problem.
Is there any faster method to achieve this interpolation? 

I want to use the EvalSome expression to solve my ode (with correct boundary conditions according to the given domain to avoid evaluation outside the domain)

Another question: is it possible to append some code as file in this forum? I think it's hard to read the code as part of the text..

-- 
You received this question notification because you are a member of
DOLFIN Team, which is an answer contact for DOLFIN.