Beispiel #1
0
void OGProjection::project_global(Hermes::Tuple<Space *> spaces, 
                                  Hermes::Tuple< std::pair<WeakForm::matrix_form_val_t, 
                                  WeakForm::matrix_form_ord_t> > proj_biforms, 
                                  Hermes::Tuple< std::pair<WeakForm::vector_form_val_t, 
                                  WeakForm::vector_form_ord_t> > proj_liforms, 
                                  Hermes::Tuple<MeshFunction*> source_meshfns, 
                                  scalar* target_vec, MatrixSolverType matrix_solver)
{
  _F_
  unsigned int n = spaces.size();
  unsigned int n_biforms = proj_biforms.size();
  if (n_biforms == 0)
    error("Please use the simpler version of project_global with the argument Hermes::Tuple<ProjNormType> proj_norms if you do not provide your own projection norm.");
  if (n_biforms != proj_liforms.size())
    error("Mismatched numbers of projection forms in project_global().");
  if (n != n_biforms)
    error("Mismatched numbers of projected functions and projection forms in project_global().");

  // This is needed since spaces may have their DOFs enumerated only locally
  // when they come here.
  int ndof = Space::assign_dofs(spaces);

  // Define projection weak form.
  WeakForm* proj_wf = new WeakForm(n);
  for (unsigned int i = 0; i < n; i++) {
    proj_wf->add_matrix_form(i, i, proj_biforms[i].first, proj_biforms[i].second);
    proj_wf->add_vector_form(i, proj_liforms[i].first, proj_liforms[i].second,
                    HERMES_ANY, source_meshfns[i]);
  }

  project_internal(spaces, proj_wf, target_vec, matrix_solver);
}
Beispiel #2
0
// Source function.
void source_fn(int n, Hermes::Tuple<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]);
  }
}
Beispiel #3
0
// Set Dirichlet boundary values.
void ModuleBasic::set_dirichlet_values(const std::vector<int> &bdy_markers_dirichlet,
                                       const std::vector<double> &bdy_values_dirichlet)
{
  Hermes::Tuple<int> tm;
  tm = bdy_markers_dirichlet;
  Hermes::Tuple<double> tv;
  tv = bdy_values_dirichlet;
  if (tm.size() != tv.size()) error("Mismatched numbers of Dirichlet boundary markers and values.");
  for (unsigned int i = 0; i < tm.size(); i++) this->bc_values.add_const(tm[i], tv[i]);
}
Beispiel #4
0
void omega_dt_fn(int n, Hermes::Tuple<scalar*> values, Hermes::Tuple<scalar*> dx, Hermes::Tuple<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
  }
}
Beispiel #5
0
Adapt::Adapt(Hermes::Tuple< Space* > spaces_, Hermes::Tuple<ProjNormType>
        proj_norms) :
    num_act_elems(-1),
    have_errors(false),
    have_coarse_solutions(false),
    have_reference_solutions(false)
{
  // sanity check
  if (proj_norms.size() > 0 && spaces_.size() != proj_norms.size()) 
    error("Mismatched numbers of spaces and projection types in Adapt::Adapt().");

  this->num = spaces_.size();

  // sanity checks
  error_if(this->num <= 0, "Too few components (%d), only %d supported.", this->num, H2D_MAX_COMPONENTS);
  error_if(this->num >= H2D_MAX_COMPONENTS, "Too many components (%d), only %d supported.", this->num, H2D_MAX_COMPONENTS);
  for (int i = 0; i < this->num; i++) {
    if (spaces_[i] == NULL) error("spaces[%d] is NULL in Adapt::Adapt().", i);
    this->spaces.push_back(spaces_[i]); 
  }

  // reset values
  memset(errors, 0, sizeof(errors));
  memset(form, 0, sizeof(form));
  memset(ord, 0, sizeof(ord));
  memset(sln, 0, sizeof(sln));
  memset(rsln, 0, sizeof(rsln));

  if (proj_norms.size() > 0) {
    for (int i = 0; i < this->num; i++) {
      switch (proj_norms[i]) {
        case HERMES_L2_NORM: form[i][i] = l2_form<double, scalar>; ord[i][i]  = l2_form<Ord, Ord>; break;
        case HERMES_H1_NORM: form[i][i] = h1_form<double, scalar>; ord[i][i]  = h1_form<Ord, Ord>; break;
        case HERMES_H1_SEMINORM: form[i][i] = h1_semi_form<double, scalar>; ord[i][i]  = h1_semi_form<Ord, Ord>; break;
        case HERMES_HCURL_NORM: form[i][i] = hcurl_form<double, scalar>; ord[i][i]  = hcurl_form<Ord, Ord>; break;
        case HERMES_HDIV_NORM: form[i][i] = hdiv_form<double, scalar>; ord[i][i]  = hdiv_form<Ord, Ord>; break;
        default: error("Unknown projection type in Adapt::Adapt().");
      }
    }
  }
}
Beispiel #6
0
double error_total(double (*efn)(MeshFunction*, MeshFunction*, RefMap*, RefMap*),
       double (*nfn)(MeshFunction*, RefMap*), Hermes::Tuple<Solution*>& slns1, Hermes::Tuple<Solution*>& slns2  )
{
  double error = 0.0, norm = 0.0;

  for (unsigned int i=0; i < slns1.size(); i++) {
    error += sqr(calc_abs_error(efn, slns1[i], slns2[i]));
    if (nfn) norm += sqr(calc_norm(nfn, slns2[i]));
  }

  return (nfn ? sqrt(error/norm) : sqrt(error));
}
Beispiel #7
0
void OGProjection::project_internal(Hermes::Tuple<Space *> spaces, WeakForm* wf, 
                                    scalar* target_vec, MatrixSolverType matrix_solver)
{
  _F_
  unsigned int n = spaces.size();

  // sanity checks
  if (n <= 0 || n > 10) error("Wrong number of projected functions in project_internal().");
  for (unsigned int i = 0; i < n; i++) if(spaces[i] == NULL) error("this->spaces[%d] == NULL in project_internal().", i);
  if (spaces.size() != n) error("Number of spaces must match number of projected functions in project_internal().");

  // this is needed since spaces may have their DOFs enumerated only locally.
  int ndof = Space::assign_dofs(spaces);

  // Initialize DiscreteProblem.
  bool is_linear = true;
  DiscreteProblem* dp = new DiscreteProblem(wf, spaces, is_linear);

  SparseMatrix* matrix = create_matrix(matrix_solver);
  Vector* rhs = create_vector(matrix_solver);
  Solver* solver = create_linear_solver(matrix_solver, matrix, rhs);

  dp->assemble(matrix, rhs, false);

  // Calculate the coefficient vector.
  bool solved = solver->solve();
  scalar* coeffs;
  if (solved) 
    coeffs = solver->get_solution();

  if (target_vec != NULL) 
    for (int i=0; i<ndof; i++) target_vec[i] = coeffs[i];
    
  delete solver;
  delete matrix;
  delete rhs;
  delete dp;
  delete wf;
}
Beispiel #8
0
bool calc_errors(Hermes::Tuple<Solution* > left, Hermes::Tuple<Solution *> right, Hermes::Tuple<double> & err_abs, Hermes::Tuple<double> & norm_vals, 
                 double & err_abs_total, double & norm_total, double & err_rel_total, Hermes::Tuple<ProjNormType> norms)
{
  bool default_norms = false;
  // Checks.
  if(left.size() != right.size())
    return false;
  if (norms != Hermes::Tuple<ProjNormType>())
  {
    if(left.size() != norms.size())
      return false;
  }
  else
    default_norms = true;
  
  // Zero the resulting Tuples.
  err_abs.clear();
  norm_vals.clear();

  // Zero the sums.
  err_abs_total = 0;
  norm_total = 0;
  err_rel_total = 0;
  
  // Calculation.
  for(unsigned int i = 0; i < left.size(); i++)
  {
    err_abs.push_back(calc_abs_error(left[i], right[i], default_norms ? HERMES_H1_NORM : norms[i]));
    norm_vals.push_back(calc_norm(right[i], default_norms ? HERMES_H1_NORM : norms[i]));
    err_abs_total += err_abs[i] * err_abs[i];
    norm_total += norm_vals[i] * norm_vals[i];
  }

  err_abs_total = sqrt(err_abs_total);
  norm_total = sqrt(norm_total);
  err_rel_total = err_abs_total / norm_total * 100.;
  
  // Everything went well, return appropriate flag.
  return true;
}
Beispiel #9
0
void OGProjection::project_global(Hermes::Tuple<Space *> spaces, Hermes::Tuple<Solution *> sols_src, 
                                  Hermes::Tuple<Solution *> sols_dest, MatrixSolverType matrix_solver, 
                                  Hermes::Tuple<ProjNormType> proj_norms)
{
  _F_
  
  scalar* target_vec = new scalar[Space::get_num_dofs(spaces)];
  Hermes::Tuple<MeshFunction *> ref_slns_mf;
  for (unsigned int i = 0; i < sols_src.size(); i++) 
    ref_slns_mf.push_back(static_cast<MeshFunction*>(sols_src[i]));
  
  OGProjection::project_global(spaces, ref_slns_mf, target_vec, matrix_solver, proj_norms);

  Solution::vector_to_solutions(target_vec, spaces, sols_dest);
  
  delete [] target_vec;
}
Beispiel #10
0
// Filter for entropy which uses the constants defined above.
static void calc_entropy_estimate_func(int n, Hermes::Tuple<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));
};
Beispiel #11
0
int main(int argc, char* argv[])
{
  if (NUMBER_OF_EIGENVALUES > 6) error("Maximum number of eigenvalues is 6.");
  info("Desired number of eigenvalues: %d.", NUMBER_OF_EIGENVALUES);

  // Load the mesh.
  Mesh mesh;
  H2DReader mloader;
  mloader.load("domain.mesh", &mesh);

  // Perform initial mesh refinements (optional).
  for (int i = 0; i < INIT_REF_NUM; i++) mesh.refine_all_elements();

  // Enter boundary markers. 
  // Note: "essential" means that solution value is prescribed.
  BCTypes bc_types;
  bc_types.add_bc_dirichlet(Hermes::Tuple<int>(BDY_BOTTOM, BDY_RIGHT, BDY_TOP, BDY_LEFT));

  // Enter Dirichlet boudnary values.
  BCValues bc_values;
  bc_values.add_zero(Hermes::Tuple<int>(BDY_BOTTOM, BDY_RIGHT, BDY_TOP, BDY_LEFT));

  // Create an H1 space with default shapeset.
  H1Space space(&mesh, &bc_types, &bc_values, P_INIT);

  // Initialize the weak formulation for the left hand side i.e. H 
  WeakForm wf_left, wf_right;
  wf_left.add_matrix_form(callback(bilinear_form_left));
  wf_right.add_matrix_form(callback(bilinear_form_right));

  // Initialize refinement selector.
  H1ProjBasedSelector selector(CAND_LIST, CONV_EXP, H2DRS_DEFAULT_ORDER);

  // Initialize views.
  ScalarView sview_1("", new WinGeom(0, 0, 350, 250));
  sview_1.show_mesh(false);
  sview_1.fix_scale_width(60);
  ScalarView sview_2("", new WinGeom(360, 0, 350, 250));
  sview_2.show_mesh(false);
  sview_2.fix_scale_width(60);
  ScalarView sview_3("", new WinGeom(720, 0, 350, 250));
  sview_3.show_mesh(false);
  sview_3.fix_scale_width(60);
  ScalarView sview_4("", new WinGeom(0, 305, 350, 250));
  sview_4.show_mesh(false);
  sview_4.fix_scale_width(60);
  ScalarView sview_5("", new WinGeom(360, 305, 350, 250));
  sview_5.show_mesh(false);
  sview_5.fix_scale_width(60);
  ScalarView sview_6("", new WinGeom(720, 305, 350, 250));
  sview_6.show_mesh(false);
  sview_6.fix_scale_width(60);
  OrderView  oview("Polynomial orders", new WinGeom(1080, 0, 410, 350));

  // DOF and CPU convergence graphs.
  SimpleGraph graph_dof_est, graph_cpu_est;

  // Time measurement.
  TimePeriod cpu_time;
  cpu_time.tick();

  Solution sln[NUMBER_OF_EIGENVALUES], ref_sln[NUMBER_OF_EIGENVALUES];

  // Adaptivity loop:
  int as = 1;
  bool done = false;
  do
  {
    info("---- Adaptivity step %d:", as);
    info("Solving on reference mesh.");

    // Construct globally refined reference mesh and setup reference space.
    Space* ref_space = construct_refined_space(&space);
    int ref_ndof = Space::get_num_dofs(ref_space);
    info("ref_ndof: %d.", ref_ndof);

    // Initialize matrices and matrix solver on referenc emesh.
    SparseMatrix* matrix_left = create_matrix(matrix_solver);
    SparseMatrix* matrix_right = create_matrix(matrix_solver);
    Solver* solver = create_linear_solver(matrix_solver, matrix_left);

    // Assemble the matrices on reference mesh.
    bool is_linear = true;
    DiscreteProblem* dp_left = new DiscreteProblem(&wf_left, ref_space, is_linear);
    dp_left->assemble(matrix_left);
    DiscreteProblem* dp_right = new DiscreteProblem(&wf_right, ref_space, is_linear);
    dp_right->assemble(matrix_right);

    // Time measurement.
    cpu_time.tick();

    // Write matrix_left in MatrixMarket format.
    write_matrix_mm("mat_left.mtx", matrix_left);

    // Write matrix_left in MatrixMarket format.
    write_matrix_mm("mat_right.mtx", matrix_right);

    // Time measurement.
    cpu_time.tick(HERMES_SKIP);

    // Calling Python eigensolver. Solution will be written to "eivecs.dat".
    info("Calling Pysparse...");
    char call_cmd[255];
    sprintf(call_cmd, "python solveGenEigenFromMtx.py mat_left.mtx mat_right.mtx %g %d %g %d", 
	    TARGET_VALUE, NUMBER_OF_EIGENVALUES, TOL, MAX_ITER);
    system(call_cmd);
    info("Pysparse finished.");

    // Initializing solution vector, solution and ScalarView.
    double* ref_coeff_vec = new double[ref_ndof];
    //Solution sln[NUMBER_OF_EIGENVALUES], ref_sln[NUMBER_OF_EIGENVALUES];
    //ScalarView view("Solution", new WinGeom(0, 0, 440, 350));

    // Reading solution vectors from file and visualizing.
    double eigenval[NUMBER_OF_EIGENVALUES];
    FILE *file = fopen("eivecs.dat", "r");
    char line [64];                  // Maximum line size.
    fgets(line, sizeof line, file);  // ref_ndof
    int n = atoi(line);            
    if (n != ref_ndof) error("Mismatched ndof in the eigensolver output file.");  
    fgets(line, sizeof line, file);  // Number of eigenvectors in the file.
    int neig = atoi(line);
    if (neig != NUMBER_OF_EIGENVALUES) error("Mismatched number of eigenvectors in the eigensolver output file.");  
    for (int ieig = 0; ieig < NUMBER_OF_EIGENVALUES; ieig++) {
      // Get next eigenvalue from the file
      fgets(line, sizeof line, file);  // eigenval
      eigenval[ieig] = atof(line);            
      // Get the corresponding eigenvector.
      for (int i = 0; i < ref_ndof; i++) {  
        fgets(line, sizeof line, file);
        ref_coeff_vec[i] = atof(line);
      }

      // Convert coefficient vector into a Solution.
      Solution::vector_to_solution(ref_coeff_vec, ref_space, &(ref_sln[ieig]));

      // Project the fine mesh solution onto the coarse mesh.
      info("Projecting reference solution %d on coarse mesh.", ieig);
      OGProjection::project_global(&space, &(ref_sln[ieig]), &(sln[ieig]), matrix_solver);
    }  
    fclose(file);
    delete [] ref_coeff_vec;

    // FIXME: Below, the adaptivity is done for the last eigenvector only,
    // this needs to be changed to take into account all eigenvectors.

    // View the coarse mesh solution and polynomial orders.
    
    char title[100];
    if (NUMBER_OF_EIGENVALUES > 0) {
      sprintf(title, "Solution 0, val = %g", eigenval[0]);
      sview_1.set_title(title);
      sview_1.show(&(sln[0]));
    }
    if (NUMBER_OF_EIGENVALUES > 1) {
      sprintf(title, "Solution 1, val = %g", eigenval[1]);
      sview_2.set_title(title);
      sview_2.show(&(sln[1]));
    }
    if (NUMBER_OF_EIGENVALUES > 2) {
      sprintf(title, "Solution 2, val = %g", eigenval[2]);
      sview_3.set_title(title);
      sview_3.show(&(sln[2]));
    }
    if (NUMBER_OF_EIGENVALUES > 3) {
      sprintf(title, "Solution 3, val = %g", eigenval[3]);
      sview_4.set_title(title);
      sview_4.show(&(sln[3]));
    }
    if (NUMBER_OF_EIGENVALUES > 4) {
      sprintf(title, "Solution 4, val = %g", eigenval[4]);
      sview_5.set_title(title);
      sview_5.show(&(sln[4]));
    }
    if (NUMBER_OF_EIGENVALUES > 5) {
      sprintf(title, "Solution 5, val = %g", eigenval[5]);
      sview_6.set_title(title);
      sview_6.show(&(sln[5]));
    }
    oview.show(&space);

    // Calculate element errors and total error estimate.
    info("Calculating error estimate.");
    Hermes::Tuple<Space *> spaces;
    for(int i = 0; i < NUMBER_OF_EIGENVALUES; i++)
        spaces.push_back(&space);
    Hermes::Tuple<ProjNormType> proj_norms;
    for(int i = 0; i < NUMBER_OF_EIGENVALUES; i++)
        proj_norms.push_back(HERMES_H1_NORM);
    Adapt* adaptivity = new Adapt(spaces, proj_norms);
    bool solutions_for_adapt = true;
    
    Hermes::Tuple<Solution *> slns;
    if (NUMBER_OF_EIGENVALUES > 0) slns.push_back(&sln[0]);
    if (NUMBER_OF_EIGENVALUES > 1) slns.push_back(&sln[1]);
    if (NUMBER_OF_EIGENVALUES > 2) slns.push_back(&sln[2]);
    if (NUMBER_OF_EIGENVALUES > 3) slns.push_back(&sln[3]);
    if (NUMBER_OF_EIGENVALUES > 4) slns.push_back(&sln[4]);
    if (NUMBER_OF_EIGENVALUES > 5) slns.push_back(&sln[5]);
    
    Hermes::Tuple<Solution *> ref_slns;
    if (NUMBER_OF_EIGENVALUES > 0) ref_slns.push_back(&ref_sln[0]);
    if (NUMBER_OF_EIGENVALUES > 1) ref_slns.push_back(&ref_sln[1]);
    if (NUMBER_OF_EIGENVALUES > 2) ref_slns.push_back(&ref_sln[2]);
    if (NUMBER_OF_EIGENVALUES > 3) ref_slns.push_back(&ref_sln[3]);
    if (NUMBER_OF_EIGENVALUES > 4) ref_slns.push_back(&ref_sln[4]);
    if (NUMBER_OF_EIGENVALUES > 5) ref_slns.push_back(&ref_sln[5]);
    Hermes::Tuple<double> component_errors;
    double err_est_rel = adaptivity->calc_err_est(slns, ref_slns, solutions_for_adapt, 
                         HERMES_TOTAL_ERROR_REL | HERMES_ELEMENT_ERROR_REL, &component_errors) * 100;

    // Report results.
    info("ndof_coarse: %d, ndof_fine: %d.", Space::get_num_dofs(&space), Space::get_num_dofs(ref_space));
    if (NUMBER_OF_EIGENVALUES > 0) info("err_est_rel[0]: %g%%", component_errors[0] * 100);
    if (NUMBER_OF_EIGENVALUES > 1) info("err_est_rel[1]: %g%%", component_errors[1] * 100);
    if (NUMBER_OF_EIGENVALUES > 2) info("err_est_rel[2]: %g%%", component_errors[2] * 100);
    if (NUMBER_OF_EIGENVALUES > 3) info("err_est_rel[3]: %g%%", component_errors[3] * 100);
    if (NUMBER_OF_EIGENVALUES > 4) info("err_est_rel[4]: %g%%", component_errors[4] * 100);
    if (NUMBER_OF_EIGENVALUES > 5) info("err_est_rel[5]: %g%%", component_errors[5] * 100);
   
    // Time measurement.
    cpu_time.tick();

    // Add entry to DOF and CPU convergence graphs.
    graph_dof_est.add_values(Space::get_num_dofs(&space), err_est_rel);
    graph_dof_est.save("conv_dof_est.dat");
    graph_cpu_est.add_values(cpu_time.accumulated(), err_est_rel);
    graph_cpu_est.save("conv_cpu_est.dat");

    // If err_est too large, adapt the mesh.
    if (err_est_rel < ERR_STOP) done = true;
    else
    {
      info("Adapting coarse mesh.");
      Hermes::Tuple<RefinementSelectors::Selector *> selectors;
      for(int i = 0; i < NUMBER_OF_EIGENVALUES; i++)
        selectors.push_back(&selector);
      done = adaptivity->adapt(selectors, THRESHOLD, STRATEGY, MESH_REGULARITY);

      // Increase the counter of performed adaptivity steps.
      if (done == false)  as++;
    }
    if (Space::get_num_dofs(&space) >= NDOF_STOP) done = true;

    // Clean up.
    delete solver;
    delete matrix_left;
    delete matrix_right;
    delete adaptivity;
    if(done == false) delete ref_space->get_mesh();
    delete ref_space;
    delete dp_left;
    delete dp_right;
  }
  while (done == false);

  // Wait for all views to be closed.
  View::wait();
  return 0;
};
Beispiel #12
0
bool Adapt::adapt(Hermes::Tuple<RefinementSelectors::Selector *> refinement_selectors, double thr, int strat,
                  int regularize, double to_be_processed)
{
    error_if(!have_errors, "element errors have to be calculated first, call Adapt::calc_err_est().");
    error_if(refinement_selectors == Hermes::Tuple<RefinementSelectors::Selector *>(), "selector not provided");
    if (spaces.size() != refinement_selectors.size()) error("Wrong number of refinement selectors.");
    TimePeriod cpu_time;

    //get meshes
    int max_id = -1;
    Mesh* meshes[H2D_MAX_COMPONENTS];
    for (int j = 0; j < this->num; j++) {
        meshes[j] = this->spaces[j]->get_mesh();
        rsln[j]->set_quad_2d(&g_quad_2d_std);
        rsln[j]->enable_transform(false);
        if (meshes[j]->get_max_element_id() > max_id)
            max_id = meshes[j]->get_max_element_id();
    }

    //reset element refinement info
    AUTOLA2_OR(int, idx, max_id + 1, this->num + 1);
    for(int j = 0; j < max_id; j++)
        for(int l = 0; l < this->num; l++)
            idx[j][l] = -1; // element not refined

    double err0_squared = 1000.0;
    double processed_error_squared = 0.0;

    vector<ElementToRefine> elem_inx_to_proc; //list of indices of elements that are going to be processed
    elem_inx_to_proc.reserve(num_act_elems);

    //adaptivity loop
    double error_squared_threshod = -1; //an error threshold that breaks the adaptivity loop in a case of strategy 1
    int num_exam_elem = 0; //a number of examined elements
    int num_ignored_elem = 0; //a number of ignored elements
    int num_not_changed = 0; //a number of element that were not changed
    int num_priority_elem = 0; //a number of elements that were processed using priority queue

    bool first_regular_element = true; //true if first regular element was not processed yet
    int inx_regular_element = 0;
    while (inx_regular_element < num_act_elems || !priority_queue.empty())
    {
        int id, comp, inx_element;

        //get element identification
        if (priority_queue.empty()) {
            id = regular_queue[inx_regular_element].id;
            comp = regular_queue[inx_regular_element].comp;
            inx_element = inx_regular_element;
            inx_regular_element++;
        }
        else {
            id = priority_queue.front().id;
            comp = priority_queue.front().comp;
            inx_element = -1;
            priority_queue.pop();
            num_priority_elem++;
        }
        num_exam_elem++;

        //get info linked with the element
        double err_squared = errors[comp][id];
        Mesh* mesh = meshes[comp];
        Element* e = mesh->get_element(id);

        if (!should_ignore_element(inx_element, mesh, e)) {
            //check if adaptivity loop should end
            if (inx_element >= 0) {
                //prepare error threshold for strategy 1
                if (first_regular_element) {
                    error_squared_threshod = thr * err_squared;
                    first_regular_element = false;
                }

                // first refinement strategy:
                // refine elements until prescribed amount of error is processed
                // if more elements have similar error refine all to keep the mesh symmetric
                if ((strat == 0) && (processed_error_squared > sqrt(thr) * errors_squared_sum)
                        && fabs((err_squared - err0_squared)/err0_squared) > 1e-3) break;

                // second refinement strategy:
                // refine all elements whose error is bigger than some portion of maximal error
                if ((strat == 1) && (err_squared < error_squared_threshod)) break;

                if ((strat == 2) && (err_squared < thr)) break;

                if ((strat == 3) &&
                        ( (err_squared < error_squared_threshod) ||
                          ( processed_error_squared > 1.5 * to_be_processed )) ) break;
            }

            // get refinement suggestion
            ElementToRefine elem_ref(id, comp);
            int current = this->spaces[comp]->get_element_order(id);
            bool refined = refinement_selectors[comp]->select_refinement(e, current, rsln[comp], elem_ref);

            //add to a list of elements that are going to be refined
            if (can_refine_element(mesh, e, refined, elem_ref) ) {
                idx[id][comp] = (int)elem_inx_to_proc.size();
                elem_inx_to_proc.push_back(elem_ref);
                err0_squared = err_squared;
                processed_error_squared += err_squared;
            }
            else {
                debug_log("Element (id:%d, comp:%d) not changed", e->id, comp);
                num_not_changed++;
            }
        }
        else {
            num_ignored_elem++;
        }
    }

    verbose("Examined elements: %d", num_exam_elem);
    verbose(" Elements taken from priority queue: %d", num_priority_elem);
    verbose(" Ignored elements: %d", num_ignored_elem);
    verbose(" Not changed elements: %d", num_not_changed);
    verbose(" Elements to process: %d", elem_inx_to_proc.size());
    bool done = false;
    if (num_exam_elem == 0)
        done = true;
    else if (elem_inx_to_proc.empty())
    {
        warn("None of the elements selected for refinement could be refined. Adaptivity step not successful, returning 'true'.");
        done = true;
    }

    //fix refinement if multimesh is used
    fix_shared_mesh_refinements(meshes, elem_inx_to_proc, idx, refinement_selectors);

    //apply refinements
    apply_refinements(elem_inx_to_proc);

    // in singlemesh case, impose same orders across meshes
    homogenize_shared_mesh_orders(meshes);

    // mesh regularization
    if (regularize >= 0)
    {
        if (regularize == 0)
        {
            regularize = 1;
            warn("Total mesh regularization is not supported in adaptivity. 1-irregular mesh is used instead.");
        }
        for (int i = 0; i < this->num; i++)
        {
            int* parents;
            parents = meshes[i]->regularize(regularize);
            this->spaces[i]->distribute_orders(meshes[i], parents);
            delete [] parents;
        }
    }

    for (int j = 0; j < this->num; j++)
        rsln[j]->enable_transform(true);

    verbose("Refined elements: %d", elem_inx_to_proc.size());
    report_time("Refined elements in: %g s", cpu_time.tick().last());

    //store for the user to retrieve
    last_refinements.swap(elem_inx_to_proc);

    have_errors = false;
    if (strat == 2 && done == true)
        have_errors = true; // space without changes

    // since space changed, assign dofs:
    Space::assign_dofs(this->spaces);

    return done;
}
Beispiel #13
0
void OGProjection::project_global(Hermes::Tuple<Space *> spaces, Hermes::Tuple<MeshFunction*> source_meshfns, 
                   scalar* target_vec, MatrixSolverType matrix_solver, Hermes::Tuple<ProjNormType> proj_norms)
{
  _F_
  int n = spaces.size();  

  // define temporary projection weak form
  WeakForm* proj_wf = new WeakForm(n);
  int found[100];
  for (int i = 0; i < 100; i++) found[i] = 0;
  for (int i = 0; i < n; i++) 
  {
    int norm;
    if (proj_norms == Hermes::Tuple<ProjNormType>()) {
      ESpaceType space_type = spaces[i]->get_type();
      switch (space_type) {
        case HERMES_H1_SPACE: norm = HERMES_H1_NORM; break;
        case HERMES_HCURL_SPACE: norm = HERMES_HCURL_NORM; break;
        case HERMES_HDIV_SPACE: norm = HERMES_HDIV_NORM; break;
        case HERMES_L2_SPACE: norm = HERMES_L2_NORM; break;
        default: error("Unknown space type in OGProjection::project_global().");
      }
    }
    else norm = proj_norms[i];
    if (norm == HERMES_H1_NORM) 
    {
      found[i] = 1;
      proj_wf->add_matrix_form(i, i, H1projection_biform<double, scalar>, H1projection_biform<Ord, Ord>);
      proj_wf->add_vector_form(i, H1projection_liform<double, scalar>, H1projection_liform<Ord, Ord>,
                               HERMES_ANY, source_meshfns[i]);
    }
    if (norm == HERMES_H1_SEMINORM) 
    {
      found[i] = 1;
      proj_wf->add_matrix_form(i, i, H1_semi_projection_biform<double, scalar>, H1_semi_projection_biform<Ord, Ord>);
      proj_wf->add_vector_form(i, H1_semi_projection_liform<double, scalar>, H1_semi_projection_liform<Ord, Ord>,
                               HERMES_ANY, source_meshfns[i]);
    }
    if (norm == HERMES_HCURL_NORM) 
    {
      found[i] = 1;
      proj_wf->add_matrix_form(i, i, Hcurlprojection_biform<double, scalar>, Hcurlprojection_biform<Ord, Ord>);
      proj_wf->add_vector_form(i, Hcurlprojection_liform<double, scalar>, Hcurlprojection_liform<Ord, Ord>,
                               HERMES_ANY, source_meshfns[i]);
    }
    if (norm == HERMES_HDIV_NORM) 
    {
      found[i] = 1;
      proj_wf->add_matrix_form(i, i, Hdivprojection_biform<double, scalar>, Hdivprojection_biform<Ord, Ord>);
      proj_wf->add_vector_form(i, Hdivprojection_liform<double, scalar>, Hdivprojection_liform<Ord, Ord>,
                               HERMES_ANY, source_meshfns[i]);
    }
    if (norm == HERMES_L2_NORM) 
    {
      found[i] = 1;
      proj_wf->add_matrix_form(i, i, L2projection_biform<double, scalar>, L2projection_biform<Ord, Ord>);
      proj_wf->add_vector_form(i, L2projection_liform<double, scalar>, L2projection_liform<Ord, Ord>,
                               HERMES_ANY, source_meshfns[i]);
    }
  }
  for (int i=0; i < n; i++) 
  {
    if (found[i] == 0) 
    {
      warn("index of component: %d\n", i);
      error("Wrong projection norm in project_global().");
    }
  }

  project_internal(spaces, proj_wf, target_vec, matrix_solver);
}