On 12/21/2017 05:54 AM, Tongzhao Gong wrote:
I'm using the deal.ii library to solve a phase-field model to simulate
columnar dendrite growth during solidification, and several references
proposed a "moving frame method" to improve the computational efficiency. For
example, I employed a 100*500 rectangle computational mesh and the dendrite
will growth in the y direction. Once the distance between the dendrite tip and
the top boundary is less than 100, the bottom region with a height of 100 is
shifted and cut off. The value of the unknown fields in the shifted region are
frozen and saved, while a new region with a height of 100 should be added onto
the top boundary, and the fields in this new region at this moving time step
is konwn. If this "moving frame method" could be realized by the deal.II library?
Yes, though it will be easier to think of it as "moving the solution down by
100 cells", rather than removing 100 cells at the bottom and adding 100 cells
at the top. In a first step, it may also be easiest to work with uniform meshes.
The way I would implement this is as follows: When you move the solution down
by 100 cells (say, by a distance Y corresponding to 100 cells), then that
means interpolating the old solution onto the mesh but shifted. You would do
this using VectorTools::interpolate with a function object of the following
form (pseudo-code, you'll need to fill in the details):
class ShiftedSolution : Function<dim> {
double value (const Point<dim> &p, ...) const {
if (p[1]+Y < top boundary)
return solution at (p[0], p[1]+Y);
else
return whatever value you assume exists beyond the top boundary;
}
};
For this, you need to be able to evaluate the previous solution at arbitrary
points (p[0], p[1]+Y). This is easy, however: you can use the FEFieldFunction
class for this.
This interpolation step will not be particularly cheap. But since you don't do
it very often, that probably doesn't matter much.
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
For more options, visit https://groups.google.com/d/optout.