Regarding a change in step-26 initial condition

35 views
Skip to first unread message

Pawan Kumar

unread,
Mar 13, 2020, 8:08:45 AM3/13/20
to deal.II User Group
Dear all;

To incorporate adaptive mesh refinement after each time steps in my ongoing work, I am trying to make follwoing small changes in step-26:

i. Initial condition 
ii. A rectangular domain
iii. Homogeneous Neumann BC.

But I am getting some errors(negative values) in the obtained initial condition plot. The changes I have done in the following ways(also the code and initial plot is attached for your reference):

The initial values class is defined and called similar to step 23/31:

Template <int dim>
class InitialValues : public Function<dim>
{
public:
	InitialValues()
: Function<dim>(1)
  {}
	virtual double value(const Point<dim> & p,
			const unsigned int component = 0) const;
};
template <int dim>
double InitialValues<dim>::value(const Point<dim> &p,
		const unsigned int component) const
		{
	Assert (component == 0, ExcInternalError());
	const double radius = 0.2;
	const double sigma = 0.1;

	const double centerx = 0.0;
	const double centery = 0.0;
	const double dx =  pow(p[0]-centerx,2);
	const double dy =  pow(p[1]-centery,2);
	const double distance = sqrt(dx+dy);

	double f_val = exp(-(dx+dy)/(2*pow(sigma,2)));

	if ((distance-radius)<1e-25)
		return f_val;
	else
		return 0;

		}
VectorTools::project(dof_handler,
			constraints,
			QGauss<dim>(fe.degree + 1),
			InitialValues<dim>(),
			old_solution);


The grid is generation:

GridGenerator::hyper_cube(triangulation,-1,1);

Due to homogeneous Neumann BC, I have removed the boundary class and related functions.

I could not understand the reason for this, could someone please guide me in this regard?

Thanks & Regards

Pawan



step-26.cc
visit0001.png

Wolfgang Bangerth

unread,
Mar 13, 2020, 10:55:34 AM3/13/20
to dea...@googlegroups.com
On 3/13/20 6:08 AM, Pawan Kumar wrote:
>
> To incorporate adaptive mesh refinement after each time steps in my ongoing work, I am trying to make follwoing small changes in step-26:
>
>
> i. Initial condition
>
> ii. A rectangular domain
>
> iii. Homogeneous Neumann BC.
>
>
> But I am getting some errors(negative values) in the obtained initial condition plot.

That is to be expected. The *projection* of a non-negative function
will, in general, not be a non-negative function. The reason is
essentially the same as when you do a truncated Fourier series of a
discontinuous function: You get Gibb's phenomenon.

If these negative values bother you, try *interpolating* the initial
conditions instead of projecting them.

Best
W.

--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Pawan Kumar

unread,
Mar 13, 2020, 11:23:47 AM3/13/20
to dea...@googlegroups.com
Dear Prof. Bangerth,

Thank you very much for the explanation. Interpolation of the initial condition works pretty well!

Thanks & Regards

Pawan

--
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 dealii+un...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/255a21fb-19fa-a7c6-e8d0-17158e75c709%40colostate.edu.
Reply all
Reply to author
Forward
0 new messages