Displacement interpolation at arbitrary point

Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Thank you and regards, Kathrin

Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards, Kathrin _______________________________________________ SfePy mailing list sfepy@python.org https://mail.python.org/mm3/mailman3/lists/sfepy.python.org/

Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables
field = Field.from_args('displacement', numpy.float64, 3, omega,
approx_order=1)
u = FieldVariable('u', 'unknown', field, order=0)
v = FieldVariable('v', 'test', field, primary_var_name='u')
...
# run simulation
vec = pb.solve()
# postprocessing and saving results
output_dict = vec.create_output_dict(fill_value=None,
extend=True,
linearization=None) # u als key
# struct array displacements, strain, stress
output_dict = post_process(output_dict, pb, state, extend=True)
nodes_displacements = output_dict['u'].data
strain_tensors = output_dict['cauchy_strain'].data
stress_tensors = output_dict['cauchy_stress'].data
...
p1 = numpy.array([1120, 53.375, 60.804])
coors = numpy.array([p1])
interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
Thank you, Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin _______________________________________________ SfePy mailing list sfepy@python.org https://mail.python.org/mm3/mailman3/lists/sfepy.python.org/
SfePy mailing list sfepy@python.org https://mail.python.org/mm3/mailman3/lists/sfepy.python.org/

Hi Kathrin,
On 10/30/2017 03:59 PM, Kathrin Sobe wrote:
Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables field = Field.from_args('displacement', numpy.float64, 3, omega, approx_order=1) u = FieldVariable('u', 'unknown', field, order=0) v = FieldVariable('v', 'test', field, primary_var_name='u') ... # run simulation vec = pb.solve() # postprocessing and saving results output_dict = vec.create_output_dict(fill_value=None, extend=True, linearization=None) # u als key # struct array displacements, strain, stress output_dict = post_process(output_dict, pb, state, extend=True) nodes_displacements = output_dict['u'].data strain_tensors = output_dict['cauchy_strain'].data stress_tensors = output_dict['cauchy_stress'].data ... p1 = numpy.array([1120, 53.375, 60.804]) coors = numpy.array([p1]) interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
You could try using Problem.save_restart() (to save the state after the problem is solved) and Problem.load_restart() (to load the state instead of solving again). But this needs the full problem instance, which might be too heavy. This use case was not really considered.
r.
Thank you, Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin

Hi Robert,
thanks for the clarification.
I was trying out your idea with reloading the problem. It works, but I aggree that this is not really the solution I was looking for. Saving the whole instance is too heavy.
Probably there is another approach: For doing the interpolation outside of sfepy I would need the regression functions of some of the tetrahedral elements (basis function of the tetrahedral elements) where the evaluation points are located. Is it possible to get this in a reasonable way after the domain or problem is defined, just by giving the evaluation points?
Then I could just multiply the displacements with the basis to get the interpolation of the points. Probably these things are already used inside of evaluate_at(). I was just not able to follow everything that is done in this function. Alternatively, I can calculate the basis on my own, but I prefer to reuse the things that are already there. Problably you have good recommendation? Thank you again.
Regards, Kathrin
2017-10-31 11:44 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/30/2017 03:59 PM, Kathrin Sobe wrote:
Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables field = Field.from_args('displacement', numpy.float64, 3, omega, approx_order=1) u = FieldVariable('u', 'unknown', field, order=0) v = FieldVariable('v', 'test', field, primary_var_name='u') ... # run simulation vec = pb.solve() # postprocessing and saving results output_dict = vec.create_output_dict(fill_value=None, extend=True, linearization=None) # u als
key # struct array displacements, strain, stress output_dict = post_process(output_dict, pb, state, extend=True) nodes_displacements = output_dict['u'].data strain_tensors = output_dict['cauchy_strain'].data stress_tensors = output_dict['cauchy_stress'].data
... p1 = numpy.array([1120, 53.375, 60.804]) coors = numpy.array([p1]) interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
You could try using Problem.save_restart() (to save the state after the problem is solved) and Problem.load_restart() (to load the state instead of solving again). But this needs the full problem instance, which might be too heavy. This use case was not really considered.
r.
Thank you,
Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin

On 11/01/2017 12:04 PM, Kathrin Sobe wrote:
Hi Robert,
thanks for the clarification.
I was trying out your idea with reloading the problem. It works, but I aggree that this is not really the solution I was looking for. Saving the whole instance is too heavy.
The easiest way might be to re-create the field only, and call the evaluate_at().
Probably there is another approach: For doing the interpolation outside of sfepy I would need the regression functions of some of the tetrahedral elements (basis function of the tetrahedral elements) where the evaluation points are located. Is it possible to get this in a reasonable way after the domain or problem is defined, just by giving the evaluation points?
Then I could just multiply the displacements with the basis to get the interpolation of the points. Probably these things are already used inside of evaluate_at(). I was just not able to follow everything that is done in this function. Alternatively, I can calculate the basis on my own, but I prefer to reuse the things that are already there. Problably you have good recommendation? Thank you again.
Yes, the evaluate_at() function does two things:
- for each physical point it finds its corresponding cell (element) and a
reference element point that maps to that physical point. (get_ref_coors()) 2. for each found cell, evaluate its basis in the found reference element points, and multiply it by the DOFs in the cell nodes. (evaluate_in_rc())
So 2. is essentially what you mean, but it does not make the matrix of the evaluated basis functions explicitly - it just applies the basis to the given DOFs (like a matrix action in matrix-free solvers).
It is possible to assemble the (sparse) basis matrix explicitly, given the data from 1. and the field DOF connectivity (see self.get_econn('volume', self.region) call).
Does this help? r.
Regards, Kathrin
2017-10-31 11:44 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/30/2017 03:59 PM, Kathrin Sobe wrote:
Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables field = Field.from_args('displacement', numpy.float64, 3, omega, approx_order=1) u = FieldVariable('u', 'unknown', field, order=0) v = FieldVariable('v', 'test', field, primary_var_name='u') ... # run simulation vec = pb.solve() # postprocessing and saving results output_dict = vec.create_output_dict(fill_value=None, extend=True, linearization=None) # u als
key # struct array displacements, strain, stress output_dict = post_process(output_dict, pb, state, extend=True) nodes_displacements = output_dict['u'].data strain_tensors = output_dict['cauchy_strain'].data stress_tensors = output_dict['cauchy_stress'].data
... p1 = numpy.array([1120, 53.375, 60.804]) coors = numpy.array([p1]) interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
You could try using Problem.save_restart() (to save the state after the problem is solved) and Problem.load_restart() (to load the state instead of solving again). But this needs the full problem instance, which might be too heavy. This use case was not really considered.
r.
Thank you,
Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin

Yes, thanks. I tried to recreate the field only and it works perfect. The Method evaluate_at() requires basically only the points where I want to have the interpolation and the saved displacement results vector. Cool.
Thank you and regards, Kathrin
2017-11-01 15:21 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
On 11/01/2017 12:04 PM, Kathrin Sobe wrote:
Hi Robert,
thanks for the clarification.
I was trying out your idea with reloading the problem. It works, but I aggree that this is not really the solution I was looking for. Saving the whole instance is too heavy.
The easiest way might be to re-create the field only, and call the evaluate_at().
Probably there is another approach: For doing the interpolation outside of
sfepy I would need the regression functions of some of the tetrahedral elements (basis function of the tetrahedral elements) where the evaluation points are located. Is it possible to get this in a reasonable way after the domain or problem is defined, just by giving the evaluation points?
Then I could just multiply the displacements with the basis to get the interpolation of the points. Probably these things are already used inside of evaluate_at(). I was just not able to follow everything that is done in this function. Alternatively, I can calculate the basis on my own, but I prefer to reuse the things that are already there. Problably you have good recommendation? Thank you again.
Yes, the evaluate_at() function does two things:
- for each physical point it finds its corresponding cell (element) and a
reference element point that maps to that physical point. (get_ref_coors()) 2. for each found cell, evaluate its basis in the found reference element points, and multiply it by the DOFs in the cell nodes. (evaluate_in_rc())
So 2. is essentially what you mean, but it does not make the matrix of the evaluated basis functions explicitly - it just applies the basis to the given DOFs (like a matrix action in matrix-free solvers).
It is possible to assemble the (sparse) basis matrix explicitly, given the data from 1. and the field DOF connectivity (see self.get_econn('volume', self.region) call).
Does this help? r.
Regards,
Kathrin
2017-10-31 11:44 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/30/2017 03:59 PM, Kathrin Sobe wrote:
Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables field = Field.from_args('displacement', numpy.float64, 3, omega, approx_order=1) u = FieldVariable('u', 'unknown', field, order=0) v = FieldVariable('v', 'test', field, primary_var_name='u') ... # run simulation vec = pb.solve() # postprocessing and saving results output_dict = vec.create_output_dict(fill_value=None, extend=True, linearization=None) # u
als key # struct array displacements, strain, stress output_dict = post_process(output_dict, pb, state, extend=True) nodes_displacements = output_dict['u'].data strain_tensors = output_dict['cauchy_strain'].data stress_tensors = output_dict['cauchy_stress'].data
... p1 = numpy.array([1120, 53.375, 60.804]) coors = numpy.array([p1]) interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
You could try using Problem.save_restart() (to save the state after the problem is solved) and Problem.load_restart() (to load the state instead of solving again). But this needs the full problem instance, which might be too heavy. This use case was not really considered.
r.
Thank you,
Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
I was looking for a postprocessing method to calculate the displacements for points on the surface of the mesh that are not necessarily node points. Precisely I need a interpolation function that gives me the displacement of these points. The input would be the already calculated linear elastic problem from my earlier example with the node displacements as result.
In the sfepy documentation I found the probing described here: http://sfepy.org/doc-devel/primer.html#sec-primer
Is this the right way to solve the interpolation? Do you recommend something different?
Yes, the best is to use FieldVariable.evaluate_at(), which can
evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin

On 11/01/2017 05:20 PM, Kathrin Sobe wrote:
Yes, thanks. I tried to recreate the field only and it works perfect. The Method evaluate_at() requires basically only the points where I want to have the interpolation and the saved displacement results vector. Cool.
OK, hth!
r.
Thank you and regards, Kathrin
2017-11-01 15:21 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
On 11/01/2017 12:04 PM, Kathrin Sobe wrote:
Hi Robert,
thanks for the clarification.
I was trying out your idea with reloading the problem. It works, but I aggree that this is not really the solution I was looking for. Saving the whole instance is too heavy.
The easiest way might be to re-create the field only, and call the evaluate_at().
Probably there is another approach: For doing the interpolation outside of
sfepy I would need the regression functions of some of the tetrahedral elements (basis function of the tetrahedral elements) where the evaluation points are located. Is it possible to get this in a reasonable way after the domain or problem is defined, just by giving the evaluation points?
Then I could just multiply the displacements with the basis to get the interpolation of the points. Probably these things are already used inside of evaluate_at(). I was just not able to follow everything that is done in this function. Alternatively, I can calculate the basis on my own, but I prefer to reuse the things that are already there. Problably you have good recommendation? Thank you again.
Yes, the evaluate_at() function does two things:
- for each physical point it finds its corresponding cell (element) and a
reference element point that maps to that physical point. (get_ref_coors()) 2. for each found cell, evaluate its basis in the found reference element points, and multiply it by the DOFs in the cell nodes. (evaluate_in_rc())
So 2. is essentially what you mean, but it does not make the matrix of the evaluated basis functions explicitly - it just applies the basis to the given DOFs (like a matrix action in matrix-free solvers).
It is possible to assemble the (sparse) basis matrix explicitly, given the data from 1. and the field DOF connectivity (see self.get_econn('volume', self.region) call).
Does this help? r.
Regards,
Kathrin
2017-10-31 11:44 GMT+01:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/30/2017 03:59 PM, Kathrin Sobe wrote:
Hi Robert,
thank you for the help.
The interplation was working according to the examples from tests.
This is a part of the code from my linear elastic problem, the result extraction (displacements, stress and strain) and the additional interpolation of certain points at the end:
# field and field variables field = Field.from_args('displacement', numpy.float64, 3, omega, approx_order=1) u = FieldVariable('u', 'unknown', field, order=0) v = FieldVariable('v', 'test', field, primary_var_name='u') ... # run simulation vec = pb.solve() # postprocessing and saving results output_dict = vec.create_output_dict(fill_value=None, extend=True, linearization=None) # u
als key # struct array displacements, strain, stress output_dict = post_process(output_dict, pb, state, extend=True) nodes_displacements = output_dict['u'].data strain_tensors = output_dict['cauchy_strain'].data stress_tensors = output_dict['cauchy_stress'].data
... p1 = numpy.array([1120, 53.375, 60.804]) coors = numpy.array([p1]) interpolations = u.evaluate_at(coors)
My question is: Is it possible to run the point interpolation after everything is finished and closed just by loading the saved result data. What I mean is, I run the FE simulation in a first step and I save the results (displacements, stress, strain). Later, I would like to load the saved data and run the point interpolation. As far as I understand the code, I would need the FieldVariable and the Field, but also the results data of the problem. Could you give some advice what would be necessary for saving to reload the required data structures for the interpolation?
You could try using Problem.save_restart() (to save the state after the problem is solved) and Problem.load_restart() (to load the state instead of solving again). But this needs the full problem instance, which might be too heavy. This use case was not really considered.
r.
Thank you,
Kathrin
2017-10-24 21:24 GMT+02:00 Robert Cimrman cimrman3@ntc.zcu.cz:
Hi Kathrin,
On 10/24/2017 06:08 PM, Kathrin Sobe wrote:
Hi Robert,
> > I was looking for a postprocessing method to calculate the > displacements > for points on the surface of the mesh that are not necessarily node > points. > Precisely I need a interpolation function that gives me the > displacement of > these points. The input would be the already calculated linear elastic > problem from my earlier example with the node displacements as result. > > In the sfepy documentation I found the probing described here: > http://sfepy.org/doc-devel/primer.html#sec-primer > > Is this the right way to solve the interpolation? Do you recommend > something different? > > > Yes, the best is to use FieldVariable.evaluate_at(), which can evaluate a FE field variable in physical points coordinates (i.e. not reference element coordinates) - check tests/test_projections.py or tests/test_mesh_interp.py. Unfortunately, there is no example showing this feature, just the tests.
r.
Thank you and regards,
Kathrin > >
participants (2)
-
Kathrin Sobe
-
Robert Cimrman