Hello all,
I have the following problem which has me going in circles. Perhaps
someone has ahieved a similar goal and can offer me some advice. The
model looks like this:
Step-1: Solve elastic problem using
fe_dim (FE_Q<dim>(1), dim),
The solution found is the displacment field.
Step-2: Solve the Poisson-like equation using
fe_1 (FE_Q<dim>(1), 1),
The solution found is the electrostatic potential.
In my case the *gradient* of the solution to Step-1 is used in the
construction of the Poisson-like equation, ie. I use the strains from
Step-1 as I assemble Step-2.
The problem is that in the Step-2 assmeble_matrix loop, by using
fe_values.get_function_gradients (solution_displacement,
solution_strain)
things get messed up because the loop over
typename DoFHandler<dim>::active_cell_iterator
cell = dof_handler.begin_active(),
endc = dof_handler.end();
for (; cell!=endc; ++cell)
if (cell->subdomain_id() == this_mpi_process)
now refers to fe_1 and not fe_dim and all sorts of errors occur (see
far, far below). The real problem here is the use of two different
finite element systems for two (uncoupled) problems.
My question is: Is there an elegant way to somehow keep a reference
to the Step-1 system solution and compute properties from it within the
assemble_system loop of Step-2 even though that refers to a different
finite element system?
This is quite an important consideration, especially if this is to be
done iteratively, say in a self-consistent model which needs to
iteratively pass solutions back-and-forth between Steps a few
hundreds/thousands times.
Brute force, ie. holding two sets of finite elements, or two
dof_handlers, or even two meshes in memory, seems to be a waste of
resources(?)
Any advice, discussion, or pointers to examples (I find none) would be
very greatly appreceated!
Best,
Toby
SOME ERRORS:
An error occurred in line <3655> of file <source/fe/fe_values.cc> in
function
void dealii::FEValues<dim, spacedim>::reinit(const typename
dealii::DoFHandler<dh_dim, spacedim>::cell_iterator&) [with int dim =
2, int spacedim = 2]
The violated condition was:
static_cast<const FiniteElementData<dim>&>(*this->fe) ==
static_cast<const FiniteElementData<dim>&>(cell->get_fe())
The name and call sequence of the exception was:
typename FEVB::ExcFEDontMatch()
Additional Information:
(none)
An error occurred in line <2599> of file <source/fe/fe_values.cc> in
function
void dealii::FEValuesBase<dim,
spacedim>::get_function_gradients(const InputVector&,
std::vector<std::vector<dealii::Tensor<1, spacedim>,
std::allocator<dealii::Tensor<1, spacedim> > >,
std::allocator<std::vector<dealii::Tensor<1, spacedim>,
std::allocator<dealii::Tensor<1, spacedim> > > > >&) const [with InputVector =
dealii::PETScWrappers::Vector, int dim = 2, int spacedim = 2]
The violated condition was:
gradients[i].size() == n_components
The name and call sequence of the exception was:
ExcDimensionMismatch(gradients[i].size(), n_components)
Additional Information:
Dimension 2 not equal to 1
--
Toby D. Young
Philosphy-Physics
Assitant Professor (Adiunkt)
Polish Academy of Sciences
Warszawa, Polska
www: http://www.ippt.gov.pl/~tyoung
skype: stenografia
_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii