Dear Wolfgang,
thank you for your swift reply! I am using Version 9.4.0, compiled from
source just last week - I checked the github as well but it seems like
there is no newer version available.
I've tried to make a minimal working example and will attach it. The code
now doesn't do anything expect for defining the class Delta_of_x that
inherits from Function, and its value() function, and then tries to
project said function into an FE vector of nodal values that I want to save
to disk. The current profile Delta(x) is in real, so I tried to change all
datatypes to double but it didn't seem to help with the linking error (and
I would want complex values for my actual problem, eventually). I'm
starting to suspect that I'm doing something wrong with the project
function, but I don't understand what if so. Either way, please find the
.cc and a makefile attached. I started from one of the examples at some
point so there is probably some unnecessary includes and similar things
left that I forgot to remove, I hope the problem is still understandable.
Thank you for your help, best
Kev
On Thursday, October 27, 2022 at 12:01:28 AM UTC+2 Wolfgang Bangerth wrote:
> On 10/26/22 12:57, Kev Se wrote:
> >
> > I've recently swapped from another FEM package to deal.ii and love it so
> > far. I recently managed so far to get a correct solution out for my
> > complex-valued, nonlinear PDE. So that's good! What I am a bit confused
> > about is the following: How do I use the VectorTools::project() type of
> > function for a general, complex-valued function?
> >
>
> Kev: It looks like the function you're trying to call is implemented but
> not instantiated for complex arguments. This wouldn't be the first
> function for which this is the case. Can you come up with as small a
> testcase that illustrates the problem? The program doesn't have to do
> anything useful other than showing that you can't link.
>
> Separately, though: Have you made sure that you are working with the
> latest version of deal.II? This may (or may not) be a bug that is
> already fixed.
>
> 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].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/6e75ddf1-4801-4ec9-ba3b-6a32c01f90afn%40googlegroups.com.
# Set the name of the project and target:
SET(TARGET "mwe")
# Declare all source files the target consists of. Here, this is only
# the one step-X.cc file, but as you expand your project you may wish
# to add other source files as well. If your project becomes much larger,
# you may want to either replace the following statement by something like
# FILE(GLOB_RECURSE TARGET_SRC "source/*.cc")
# FILE(GLOB_RECURSE TARGET_INC "include/*.h")
# SET(TARGET_SRC ${TARGET_SRC} ${TARGET_INC})
# or switch altogether to the large project CMakeLists.txt file discussed
# in the "CMake in user projects" page accessible from the "User info"
# page of the documentation.
SET(TARGET_SRC
${TARGET}.cc
)
# Usually, you will not need to modify anything beyond this point...
CMAKE_MINIMUM_REQUIRED(VERSION 3.3.0)
FIND_PACKAGE(deal.II 9.4.0 QUIET
HINTS ${deal.II_DIR} ${DEAL_II_DIR} ../ ../../ $ENV{DEAL_II_DIR}
)
IF(NOT ${deal.II_FOUND})
MESSAGE(FATAL_ERROR "\n"
"*** Could not locate a (sufficiently recent) version of deal.II. ***\n\n"
"You may want to either pass a flag -DDEAL_II_DIR=/path/to/deal.II to
cmake\n"
"or set an environment variable \"DEAL_II_DIR\" that contains this path."
)
ENDIF()
DEAL_II_INITIALIZE_CACHED_VARIABLES()
PROJECT(${TARGET})
DEAL_II_INVOKE_AUTOPILOT()
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/function.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/solver_gmres.h>
#include <deal.II/lac/precondition.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_out.h>
#include <deal.II/grid/grid_refinement.h>
#include <deal.II/lac/sparse_direct.h>
#include <deal.II/fe/mapping_q1.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/numerics/vector_tools.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/numerics/data_out.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_interface_values.h>
#include <deal.II/numerics/derivative_approximation.h>
#include <deal.II/meshworker/mesh_loop.h>
#include <iostream>
#include <fstream>
typedef std::complex<double> dcomp;
const dcomp J = {0, 1};
const int degree=1;
namespace mwe
{
using namespace dealii;
template <int dim>
class Delta_of_x : public Function<dim, dcomp>
{
public:
virtual dcomp value(const Point<dim> &p,
const unsigned int component = 0) const override;
};
template <int dim>
dcomp Delta_of_x<dim>::value(const Point<dim> &p,
const unsigned int component) const
{
dcomp returner = 0.0;
if (component == 0){
returner= tanh( (2.0*M_PI*p[0] - 7.5)/2.0 );
}
return returner;
}
template <int dim>
class AdvectionProblem
{
public:
AdvectionProblem();
void run();
private:
void setup_system(const bool initial_step);
void output_results(const unsigned int cycle) const;
Triangulation<dim> triangulation;
const MappingQ1<dim> mapping;
const FE_DGQ<dim> fe;
DoFHandler<dim> dof_handler;
const QGauss<dim> quadrature;
const QGauss<dim - 1> quadrature_face;
SparsityPattern sparsity_pattern;
Vector<dcomp> Delta_profile;
};
template <int dim>
AdvectionProblem<dim>::AdvectionProblem()
: mapping()
, fe(degree)
, dof_handler(triangulation)
, quadrature(fe.tensor_degree() + 2)
, quadrature_face(fe.tensor_degree() + 2)
{}
template <int dim>
void AdvectionProblem<dim>::setup_system(const bool initial_step)
{
if (initial_step){
dof_handler.distribute_dofs(fe);
Delta_profile.reinit(dof_handler.n_dofs());
}
DynamicSparsityPattern dsp(dof_handler.n_dofs());
DoFTools::make_flux_sparsity_pattern(dof_handler, dsp);
sparsity_pattern.copy_from(dsp);
}
template <int dim>
void AdvectionProblem<dim>::output_results(const unsigned int cycle) const
{
const std::string filename = "solution-" + std::to_string(cycle) + ".vtk";
std::cout << " Writing solution to <" << filename << '>' << std::endl;
std::ofstream output(filename);
Delta_of_x<dim> DeltaX;
AffineConstraints<dcomp> constraints;
constraints.close();
/************************************* PROBLEMATIC LINE BELOW ****************************/
VectorTools::project(dof_handler, constraints, QGauss<dim>(degree + 2), DeltaX, Delta_profile);
/************************************* PROBLEMATIC LINE ABOVE ****************************/
DataOut<dim> data_out;
data_out.build_patches(mapping);
data_out.attach_dof_handler(dof_handler);
data_out.add_data_vector(Delta_profile, "Delta", DataOut<dim>::type_dof_data);
data_out.write_vtk(output);
}
template <int dim>
void AdvectionProblem<dim>::run()
{
GridGenerator::hyper_cube(triangulation, 0.0, 50/(2.0*M_PI));
triangulation.refine_global(3);
setup_system(true);
output_results(0);
}
} // namespace mwe
int main()
{
try
{
mwe::AdvectionProblem<1> dgmethod;
dgmethod.run();
}
catch (std::exception &exc)
{
std::cerr << std::endl
<< std::endl
<< "----------------------------------------------------"
<< std::endl;
std::cerr << "Exception on processing: " << std::endl
<< exc.what() << std::endl
<< "Aborting!" << std::endl
<< "----------------------------------------------------"
<< std::endl;
return 1;
}
catch (...)
{
std::cerr << std::endl
<< std::endl
<< "----------------------------------------------------"
<< std::endl;
std::cerr << "Unknown exception!" << std::endl
<< "Aborting!" << std::endl
<< "----------------------------------------------------"
<< std::endl;
return 1;
}
return 0;
}