PRISMS-PF Manual
Loading...
Searching...
No Matches
nucleus_refinement_function.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/bounding_box.h>
7#include <deal.II/base/point.h>
8#include <deal.II/grid/tria_accessor.h>
9
11#include <prismspf/core/types.h>
12
14
17
18#include <prismspf/config.h>
19
20#include <iostream>
21#include <memory>
22#include <vector>
23
25
30template <unsigned int dim>
32{
33public:
35 dealii::CellAccessor<dim>; // dealii::TriaActiveIterator<dealii::CellAccessor<dim>>;
36 using prisms::CellMarkerBase<dim>::flag;
37
43
44 bool
46 [[maybe_unused]] const SimulationTimer &time_info) const override
47 {
48 for (const Nucleus<dim> &nucleus : *nuclei_list)
49 {
51 {
52 static dealii::Point<dim> unit_corner = []()
53 {
54 dealii::Point<dim> p;
55 for (unsigned int d = 0; d < dim; ++d)
56 {
57 p[d] = 1.0;
58 }
59 return p;
60 }();
61 dealii::BoundingBox<dim> nucleus_bounding_box(
62 std::make_pair<dealii::Point<dim>, dealii::Point<dim>>(
63 dealii::Point<dim>(nucleus.location -
65 dealii::Point<dim>(nucleus.location +
67 if (cell.bounding_box().has_overlap_with(nucleus_bounding_box))
68 {
69 return true;
70 }
71 }
72 }
73 return false;
74 }
75
77 const std::vector<Nucleus<dim>> *nuclei_list;
78};
79
Base class for cell markers.
Definition cell_marker_base.h:25
The class handles the stochastic nucleation in PRISMS-PF.
Definition nucleus_refinement_function.h:32
const std::vector< Nucleus< dim > > * nuclei_list
Definition nucleus_refinement_function.h:77
NucleusRefinementFunction(const NucleationParameters &_nuc_params, const std::vector< Nucleus< dim > > &_nuclei_list)
Definition nucleus_refinement_function.h:38
bool flag(const CellIterator &cell, const SimulationTimer &time_info) const override
Definition nucleus_refinement_function.h:45
const NucleationParameters * nuc_params
Definition nucleus_refinement_function.h:76
dealii::CellAccessor< dim > CellIterator
Definition nucleus_refinement_function.h:35
Definition simulation_timer.h:13
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
Definition conditional_ostreams.cc:20
Struct that holds nucleation parameters.
Definition nucleation_parameters.h:27
bool check_active(const Nucleus< dim > &nucleus, const SimulationTimer &time_info) const
Check if a nucleus is still active based on its seed time and increment. A nucleus is considered acti...
Definition nucleation_parameters.h:169
double get_refinement_radius() const
Get the refinement radius.
Definition nucleation_parameters.h:120
This class contains mutable utilities for phase field problems.
Definition nucleus.h:23