Exemplo n.º 1
0
// Source function.
void source_fn(int n, Hermes::vector<scalar*> values, scalar* out)
{
  for (int i = 0; i < n; i++)
  {
		out[i] = (nu[1][0] * Sf[1][0] * values.at(0)[i] +
        nu[1][1] * Sf[1][1] * values.at(1)[i] +
        nu[1][2] * Sf[1][2] * values.at(2)[i] +
        nu[1][3] * Sf[1][3] * values.at(3)[i]);
  }
}
Exemplo n.º 2
0
void CustomFilterDt::filter_fn(int n, Hermes::vector<double*> values, Hermes::vector<double*> dx, Hermes::vector<double*> dy,
                        double* out, double* outdx, double* outdy)
{
  for (int i = 0; i < n; i++)
  {
    double t1 = std::max(values.at(0)[i],0.0) - 1.0;
    double t2 = t1 * beta;
    double t3 = 1.0 + t1 * alpha;
    double t4 = sqr(beta) / (2.0*Le) * exp(t2 / t3);
    double t5 = (beta / (t3 * t3));
    out[i] = t4 * t5 * values.at(1)[i];
    outdx[i] = 0.0;
    outdy[i] = 0.0; // not important
  }
}
Exemplo n.º 3
0
void omega_dt_fn(int n, Hermes::vector<scalar*> values, Hermes::vector<scalar*> dx, Hermes::vector<scalar*> dy,
                        scalar* out, scalar* outdx, scalar* outdy)
{
  for (int i = 0; i < n; i++)
  {
    scalar t1 = std::max(values.at(0)[i],0.0) - 1.0;
    scalar t2 = t1 * beta;
    scalar t3 = 1.0 + t1 * alpha;
    scalar t4 = sqr(beta) / (2.0*Le) * exp(t2 / t3);
    scalar t5 = (beta / (t3 * t3));
    out[i] = t4 * t5 * values.at(1)[i];
    outdx[i] = 0.0;
    outdy[i] = 0.0; // not important
  }
}
Exemplo n.º 4
0
 void filter_fn(int n, Hermes::vector<scalar*> values, scalar* result)
 {
   for (int i = 0; i < n; i++) {
     result[i] = 0;
     for (unsigned int j = 0; j < values.size(); j++)
       result[i] += nu[j] * Sigma_f[j] * values.at(j)[i];
   }
 } 
Exemplo n.º 5
0
/// Fission source function.
inline void source_fn(int n, Hermes::vector<scalar*> values, scalar* out)
{
  for (int i = 0; i < n; i++) {
    out[i] = 0.0;
    for_each_group(g)
      out[i] += nu[1][g] * Sf[1][g] * values.at(g)[i];
  }
}
Exemplo n.º 6
0
 void Filter<Scalar>::init(Hermes::vector<MeshFunctionSharedPtr<Scalar> > solutions)
 {
   this->num = solutions.size();
   if(num > H2D_MAX_COMPONENTS)
     throw Hermes::Exceptions::Exception("Attempt to create an instance of Filter with more than 10 MeshFunctions.");
   for(int i = 0; i < this->num; i++)
     this->sln[i] = solutions.at(i);
   this->init();
 }
Exemplo n.º 7
0
 void Filter<Scalar>::init(const Hermes::vector<MeshFunction<Scalar>*>& solutions)
 {
   this->num = solutions.size();
   if(num > 10)
     throw Hermes::Exceptions::Exception("Attempt to create an instance of Filter with more than 10 MeshFunctions.");
   for(int i = 0; i < this->num; i++)
     this->sln[i] = solutions.at(i);
   this->init();
 }
Exemplo n.º 8
0
void Filter::init(Hermes::vector<MeshFunction*> solutions)
{
	this->num = solutions.size();
	if(num > 10)
		error("Attempt to create an instance of Filter with more than 10 MeshFunctions.");
	for(int i = 0; i < this->num; i++)
		this->sln[i] = solutions.at(i);
  this->init();
}
Exemplo n.º 9
0
    SimpleFilter<Scalar>::SimpleFilter(const Hermes::vector<Solution<Scalar>*>& solutions, const Hermes::vector<int>& items)
    {
      this->num = solutions.size();
      if(this->num > 10)
        error("Attempt to create an instance of Filter with more than 10 MeshFunctions.");
      if(items.size() != (unsigned) this->num)
        if(items.size() > 0)
          error("Attempt to create an instance of SimpleFilter with different supplied number of MeshFunctions than the number of types of data used from them.");

      for(int i = 0; i < this->num; i++)
      {
        this->sln[i] = solutions.at(i);
        if(items.size() > 0)
          this->item[i] = items.at(i);
        else
          this->item[i] = H2D_FN_VAL;
      }
      this->init();
      init_components();
    }
Exemplo n.º 10
0
SimpleFilter::SimpleFilter(void (*filter_fn)(int n, Hermes::vector<scalar*> values, scalar* result),
		Hermes::vector<MeshFunction*> solutions, Hermes::vector<int> items) : filter_fn(filter_fn)
{
	this->num = solutions.size();
	if(num > 10)
		error("Attempt to create an instance of Filter with more than 10 MeshFunctions.");
	if(items.size() != (unsigned) num)
		if(items.size() > 0)
			error("Attempt to create an instance of SimpleFilter with different supplied number of MeshFunctions than the number of types of data used from them.");

	for(int i = 0; i < this->num; i++)
	{
		this->sln[i] = solutions.at(i);
		if(items.size() > 0)
			this->item[i] = items.at(i);
		else
			this->item[i] =	H2D_FN_VAL;
	}
	this->init();
	init_components();
}
Exemplo n.º 11
0
void CustomFilter::filter_fn(int n, double* x, double* y, Hermes::vector<double*> values, Hermes::vector<double*> dx, Hermes::vector<double*> dy,
                      double* out, double* outdx, double* outdy)
{
  for (int i = 0; i < n; i++)
  {
    double t1 = std::max(values.at(0)[i],0.0) - 1.0;
    double t2 = t1 * beta;
    double t3 = 1.0 + t1 * alpha;
    double t4 = sqr(beta) / (2.0*Le) * exp(t2 / t3);
    double t5 = (beta / (t3 * t3)) * values.at(1)[i];
    out[i] = t4 * values.at(1)[i];
    outdx[i] = t4 * (dx.at(1)[i] + dx.at(0)[i] * t5);
    outdy[i] = t4 * (dy.at(1)[i] + dy.at(0)[i] * t5);
  }
}
Exemplo n.º 12
0
// Filter for entropy which uses the constants defined above.
static void calc_entropy_estimate_func(int n, Hermes::vector<scalar*> scalars, scalar* result)
{
  for (int i = 0; i < n; i++)
    result[i] = std::log((calc_pressure(scalars.at(0)[i], scalars.at(1)[i], scalars.at(2)[i], scalars.at(3)[i]) / P_EXT)
    / pow((scalars.at(0)[i] / RHO_EXT), KAPPA));
};
Exemplo n.º 13
0
QList<SolutionArray *> SolutionAgros::solveSolutioArray(Hermes::vector<EssentialBCs> bcs)
{
    QTime time;

    // solution agros array
    QList<SolutionArray *> solutionArrayList;

    // load the mesh file
    mesh = readMeshFromFile(tempProblemFileName() + ".mesh");
    refineMesh(mesh, true, true);

    // create an H1 space
    Hermes::vector<Space *> space;
    // create hermes solution array
    Hermes::vector<Solution *> solution;
    // create reference solution
    Hermes::vector<Solution *> solutionReference;

    // projection norms
    Hermes::vector<ProjNormType> projNormType;

    // prepare selector
    Hermes::vector<RefinementSelectors::Selector *> selector;

    // error marker
    bool isError = false;

    RefinementSelectors::Selector *select = NULL;
    switch (adaptivityType)
    {
    case AdaptivityType_H:
        select = new RefinementSelectors::HOnlySelector();
        break;
    case AdaptivityType_P:
        select = new RefinementSelectors::H1ProjBasedSelector(RefinementSelectors::H2D_P_ANISO,
                                                              Util::config()->convExp,
                                                              H2DRS_DEFAULT_ORDER);
        break;
    case AdaptivityType_HP:
        select = new RefinementSelectors::H1ProjBasedSelector(RefinementSelectors::H2D_HP_ANISO,
                                                              Util::config()->convExp,
                                                              H2DRS_DEFAULT_ORDER);
        break;
    }

    for (int i = 0; i < numberOfSolution; i++)
    {
        space.push_back(new H1Space(mesh, &bcs[i], polynomialOrder));

        // set order by element
        for (int j = 0; j < Util::scene()->labels.count(); j++)
            if (Util::scene()->labels[j]->material != Util::scene()->materials[0])
                space.at(i)->set_uniform_order(Util::scene()->labels[j]->polynomialOrder > 0 ? Util::scene()->labels[j]->polynomialOrder : polynomialOrder,
                                               QString::number(j).toStdString());

        // solution agros array
        solution.push_back(new Solution());

        if (adaptivityType != AdaptivityType_None)
        {
            // add norm
            projNormType.push_back(Util::config()->projNormType);
            // add refinement selector
            selector.push_back(select);
            // reference solution
            solutionReference.push_back(new Solution());
        }
    }

    // check for DOFs
    if (Space::get_num_dofs(space) == 0)
    {
        m_progressItemSolve->emitMessage(QObject::tr("DOF is zero"), true);
    }
    else
    {
        for (int i = 0; i < numberOfSolution; i++)
        {
            // transient
            if (analysisType == AnalysisType_Transient)
            {
                // constant initial solution
                solution.at(i)->set_const(mesh, initialCondition);
                solutionArrayList.append(solutionArray(solution.at(i)));
            }

            // nonlinear
            if ((linearityType != LinearityType_Linear) && (analysisType != AnalysisType_Transient))
            {
                solution.at(i)->set_const(mesh, 0.0);
            }
        }

        actualTime = 0.0;

        // update time function
        Util::scene()->problemInfo()->hermes()->updateTimeFunctions(actualTime);

        m_wf->set_current_time(actualTime);
        m_wf->solution = solution;
        m_wf->delete_all();
        m_wf->registerForms();

        // emit message
        if (adaptivityType != AdaptivityType_None)
            m_progressItemSolve->emitMessage(QObject::tr("Adaptivity type: %1").arg(adaptivityTypeString(adaptivityType)), false);

        double error = 0.0;

        // solution
        int maxAdaptivitySteps = (adaptivityType == AdaptivityType_None) ? 1 : adaptivitySteps;
        int actualAdaptivitySteps = -1;
        for (int i = 0; i<maxAdaptivitySteps; i++)
        {
            // set up the solver, matrix, and rhs according to the solver selection.
            SparseMatrix *matrix = create_matrix(matrixSolver);
            Vector *rhs = create_vector(matrixSolver);
            Solver *solver = create_linear_solver(matrixSolver, matrix, rhs);

            if (adaptivityType == AdaptivityType_None)
            {
                if (analysisType != AnalysisType_Transient)
                    solve(space, solution, solver, matrix, rhs);
            }
            else
            {
                // construct globally refined reference mesh and setup reference space.
                Hermes::vector<Space *> spaceReference = *Space::construct_refined_spaces(space);

                // assemble reference problem.
                solve(spaceReference, solutionReference, solver, matrix, rhs);

                if (!isError)
                {
                    // project the fine mesh solution onto the coarse mesh.
                    OGProjection::project_global(space, solutionReference, solution, matrixSolver);

                    // Calculate element errors and total error estimate.
                    Adapt adaptivity(space, projNormType);

                    // Calculate error estimate for each solution component and the total error estimate.
                    error = adaptivity.calc_err_est(solution,
                                                    solutionReference) * 100;

                    // emit signal
                    m_progressItemSolve->emitMessage(QObject::tr("Adaptivity rel. error (step: %2/%3, DOFs: %4/%5): %1%").
                                                     arg(error, 0, 'f', 3).
                                                     arg(i + 1).
                                                     arg(maxAdaptivitySteps).
                                                     arg(Space::get_num_dofs(space)).
                                                     arg(Space::get_num_dofs(spaceReference)), false, 1);
                    // add error to the list
                    m_progressItemSolve->addAdaptivityError(error, Space::get_num_dofs(space));

                    if (error < adaptivityTolerance || Space::get_num_dofs(space) >= adaptivityMaxDOFs)
                    {
                        break;
                    }
                    if (i != maxAdaptivitySteps-1) adaptivity.adapt(selector,
                                                                    Util::config()->threshold,
                                                                    Util::config()->strategy,
                                                                    Util::config()->meshRegularity);
                    actualAdaptivitySteps = i+1;
                }

                if (m_progressItemSolve->isCanceled())
                {
                    isError = true;
                    break;
                }

                // delete reference space
                for (int i = 0; i < spaceReference.size(); i++)
                {
                    delete spaceReference.at(i)->get_mesh();
                    delete spaceReference.at(i);
                }
                spaceReference.clear();
            }

            // clean up.
            delete solver;
            delete matrix;
            delete rhs;
        }

        // delete reference solution
        for (int i = 0; i < solutionReference.size(); i++)
            delete solutionReference.at(i);
        solutionReference.clear();

        // delete selector
        if (select) delete select;
        selector.clear();

        // timesteps
        if (!isError)
        {
            SparseMatrix *matrix = NULL;
            Vector *rhs = NULL;
            Solver *solver = NULL;

            // allocate dp for transient solution
            DiscreteProblem *dpTran = NULL;
            if (analysisType == AnalysisType_Transient)
            {
                // set up the solver, matrix, and rhs according to the solver selection.
                matrix = create_matrix(matrixSolver);
                rhs = create_vector(matrixSolver);
                solver = create_linear_solver(matrixSolver, matrix, rhs);
                // solver->set_factorization_scheme(HERMES_REUSE_FACTORIZATION_COMPLETELY);

                dpTran = new DiscreteProblem(m_wf, space, true);
            }

            int timesteps = (analysisType == AnalysisType_Transient) ? floor(timeTotal/timeStep) : 1;
            for (int n = 0; n<timesteps; n++)
            {
                // set actual time
                actualTime = (n+1)*timeStep;

                // update essential bc values
                Space::update_essential_bc_values(space, actualTime);
                // update timedep values
                Util::scene()->problemInfo()->hermes()->updateTimeFunctions(actualTime);

                m_wf->set_current_time(actualTime);
                m_wf->delete_all();
                m_wf->registerForms();

                // transient
                if ((timesteps > 1) && (linearityType == LinearityType_Linear))
                    isError = !solveLinear(dpTran, space, solution,
                                           solver, matrix, rhs);
                if ((timesteps > 1) && (linearityType != LinearityType_Linear))
                    isError = !solve(space, solution,
                                     solver, matrix, rhs);

                // output
                for (int i = 0; i < numberOfSolution; i++)
                {
                    solutionArrayList.append(solutionArray(solution.at(i), space.at(i), error, actualAdaptivitySteps, (n+1)*timeStep));
                }

                if (analysisType == AnalysisType_Transient)
                    m_progressItemSolve->emitMessage(QObject::tr("Transient time step (%1/%2): %3 s").
                                                     arg(n+1).
                                                     arg(timesteps).
                                                     arg(actualTime, 0, 'e', 2), false, n+2);
                if (m_progressItemSolve->isCanceled())
                {
                    isError = true;
                    break;
                }
            }

            // clean up
            if (solver) delete solver;
            if (matrix) delete matrix;
            if (rhs) delete rhs;

            if (dpTran) delete dpTran;
        }
    }
    // delete mesh
    delete mesh;

    // delete space
    for (unsigned int i = 0; i < space.size(); i++)
    {
        // delete space.at(i)->get_mesh();
        delete space.at(i);
    }
    space.clear();

    // delete last solution
    for (unsigned int i = 0; i < solution.size(); i++)
        delete solution.at(i);
    solution.clear();

    if (isError)
    {
        for (int i = 0; i < solutionArrayList.count(); i++)
            delete solutionArrayList.at(i);
        solutionArrayList.clear();
    }
    return solutionArrayList;
}