[aspect-devel] Material models with 3-D ASCII input
Marie Kajan
marie.kajan at gmail.com
Wed Jan 31 11:39:27 PST 2018
Wolfgang,
Another option would be to make the AsciiDataInitial class a *member
> variable*, rather than a base class. Presumably you do it the way you show
> above because you want to access the elements of the AsciiDataInitial
> member functions from the evaluate() function of your material model. But
> you can also achieve this by making the AsciiDataInitial thing a member.
>
Actually, my understanding of C++ terminology may be lacking, but I think
this might be what I was trying to describe in #2 and what I tried to
implement myself in #3. I'll dump a bunch of code here and you can correct
me if I'm wrong.
The relevant parts of my material model header file:
#include <aspect/material_model/interface.h>
#include <aspect/simulator_access.h>
#include <aspect/utilities.h>
namespace aspect
{
namespace MaterialModel
{
using namespace dealii;
template <int dim>
class XXXXX: public MaterialModel::Interface<dim>, public
::aspect::SimulatorAccess<dim>
{
public:
XXXXX ();
virtual
void
initialize ();
...
private:
aspect::Utilities::AsciiDataInitial<dim> field;
...
};
}
}
And the relevant parts of my material model main code:
#include <aspect/material_model/XXXXX.h>
#include <aspect/adiabatic_conditions/interface.h>
#include <aspect/utilities.h>
#include <aspect/lateral_averaging.h>
#include <fstream>
#include <iostream>
namespace aspect
{
namespace MaterialModel
{
template <int dim>
XXXXX<dim>::XXXXX()
{}
template <int dim>
void
XXXXX<dim>::initialize()
{
std::cout << "initialize start" << std::endl;
field.initialize(1);
std::cout << "initialize end" << std::endl;
}
...and I won't go further than that because the error pops up here. This is
the AsciiDataInitial::initialize() function in utilities.cc:
template <int dim>
void
AsciiDataInitial<dim>::initialize (const unsigned int components)
{
std::cout << "ascii initialize before throw" << std::endl;
AssertThrow ((dynamic_cast<const GeometryModel::SphericalShell<dim>*>
(&this->get_geometry_model()))
|| (dynamic_cast<const GeometryModel::Chunk<dim>*>
(&this->get_geometry_model())) != 0
|| (dynamic_cast<const GeometryModel::Box<dim>*>
(&this->get_geometry_model())) != 0,
ExcMessage ("This ascii data plugin can only be used when using "
"a spherical shell, chunk or box
geometry."));
std::cout << "ascii initialize after throw" << std::endl;
...it goes on, but I commented out different parts and found that I get
this type of error for any line with '&this->' or 'this->':
[marie at lehmann aspect]$ ./aspect /home/marie/cookbooks/shell_simple_3d_test.prm
------------------------------------------------------------
-----------------
-- This is ASPECT, the Advanced Solver for Problems in Earth's ConvecTion.
-- . version 2.0.0-pre (master, 179d6da)
-- . using deal.II 8.5.0
-- . using Trilinos 12.10.1
-- . using p4est 1.1.0
-- . running in DEBUG mode
-- . running with 1 MPI process
------------------------------------------------------------
-----------------
initialize start
ascii initialize before throw
[lehmann:15947] *** Process received signal ***
[lehmann:15947] Signal: Segmentation fault (11)
[lehmann:15947] Signal code: Address not mapped (1)
[lehmann:15947] Failing at address: 0xdb0
[lehmann:15947] [ 0] /lib64/libc.so.6(+0x35250)[0x7f29709e4250]
[lehmann:15947] [ 1] ./aspect(_ZNK6aspect15SimulatorAccessILi
3EE18get_geometry_modelEv+0xf)[0xd8d0af]
[lehmann:15947] [ 2] ./aspect(_ZN6aspect9Utilities16AsciiData
InitialILi3EE10initializeEj+0x3c)[0xde85dc]
[lehmann:15947] [ 3] ./aspect(_ZN6aspect13MaterialModel19XXXX
XILi3EE10initializeEv+0xe2)[0xa39a42]
[lehmann:15947] [ 4] ./aspect(_ZN6aspect9SimulatorILi3EEC1EP1
9ompi_communicator_tRN6dealii16ParameterHandlerE+0x212b)[0xcdfebb]
[lehmann:15947] [ 5] ./aspect(main+0x7de)[0x997c0b]
[lehmann:15947] [ 6] /lib64/libc.so.6(__libc_start_
main+0xf5)[0x7f29709d0b35]
[lehmann:15947] [ 7] ./aspect[0x8464d9]
[lehmann:15947] *** End of error message ***
Segmentation fault (core dumped)
I'm not sure if this means that the AsciiDataInitial functions as written
are incompatible with this type of use, or if I could make this work by
altering something in my own code. I can dump the entire code files here if
that would make things easier to figure out.
Thanks!
Marie
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.geodynamics.org/pipermail/aspect-devel/attachments/20180131/4771c246/attachment.html>
More information about the Aspect-devel
mailing list