FIMS  v0.8.1
Loading...
Searching...
No Matches
lognormal_lpdf.hpp
Go to the documentation of this file.
1
10#ifndef LOGNORMAL_LPDF
11#define LOGNORMAL_LPDF
12
14#include "../../common/fims_vector.hpp"
15#include "../../common/def.hpp"
16
17namespace fims_distributions {
34template <typename Type>
35struct LogNormalLPDF : public DensityComponentBase<Type> {
43
47
50 virtual ~LogNormalLPDF() {}
51
63 virtual const Type evaluate() {
64 // set vector size based on input type (prior, process, or data)
65 size_t n_x = this->get_n_x();
66 // get expected value vector size
67 size_t n_expected = this->get_n_expected();
68 // setup vector for recording the log probability density function values
69 this->lpdf_vec.resize(n_x);
70 std::fill(this->lpdf_vec.begin(), this->lpdf_vec.end(),
71 static_cast<Type>(0));
72 this->lpdf = static_cast<Type>(0);
73
74 // Dimension checks
75 // TODO: fix dimension check as expected values no longer used for data
76 if (n_x != n_expected) {
77 if (n_expected == 1) {
78 n_expected = n_x;
79 } else if (n_x > n_expected) {
80 n_x = n_expected;
81 }
82 }
83
84 if (this->log_sd.size() > 1 && n_x != this->log_sd.size()) {
85 throw std::invalid_argument(
86 "LognormalLPDF::Vector index out of bounds. The size of observed "
87 "data does not equal the size of the log_sd vector. The observed "
88 "data vector is of size " +
89 fims::to_string(n_x) + " and the log_sd vector is of size " +
90 fims::to_string(this->log_sd.size()));
91 }
92
93 for (size_t i = 0; i < n_x; i++) {
94#ifdef TMB_MODEL
95 if (this->input_type == "data") {
96 // if data, check if there are any NA values and skip lpdf calculation
97 // if there are See Deroba and Miller, 2016
98 // (https://doi.org/10.1016/j.fishres.2015.12.002) for the use of
99 // lognormal constant
100 if (this->get_observed(i) != this->data_observed_values->na_value) {
101 this->lpdf_vec[i] =
102 dnorm(log(this->get_observed(i)), this->get_expected(i),
103 fims_math::exp(log_sd.get_force_scalar(i)), true) -
104 log(this->get_observed(i));
105 } else {
106 this->lpdf_vec[i] = 0;
107 }
108 } else {
109 if (this->input_type == "random_effects") {
110 // if random effects, no lognormal constant needs to be applied
111 this->lpdf_vec[i] =
112 dnorm(log(this->get_observed(i)), this->get_expected(i),
113 fims_math::exp(log_sd.get_force_scalar(i)), true);
114 } else {
115 this->lpdf_vec[i] =
116 dnorm(log(this->get_observed(i)), this->get_expected(i),
117 fims_math::exp(log_sd.get_force_scalar(i)), true) -
118 log(this->get_observed(i));
119 }
120 }
121
122 this->lpdf += this->lpdf_vec[i];
123 if (this->simulate_flag) {
124 FIMS_SIMULATE_F(this->of) { // preprocessor definition in interface.hpp
125 // this simulates data that is mean biased
126 if (this->input_type == "data") {
127 this->data_observed_values->at(i) = fims_math::exp(
128 rnorm(this->get_expected(i),
129 fims_math::exp(log_sd.get_force_scalar(i))));
130 }
131 if (this->input_type == "random_effects") {
132 (*this->re)[i] = fims_math::exp(
133 rnorm(this->get_expected(i),
134 fims_math::exp(log_sd.get_force_scalar(i))));
135 }
136 if (this->input_type == "prior") {
137 (*(this->priors[i]))[0] = fims_math::exp(
138 rnorm(this->get_expected(i),
139 fims_math::exp(log_sd.get_force_scalar(i))));
140 }
141 }
142 }
143#endif
144 }
145#ifdef TMB_MODEL
146 vector<Type> lognormal_observed_values = this->observed_values.to_tmb();
147 // FIMS_REPORT_F(lognormal_observed_values, this->of);
148#endif
149 return (this->lpdf);
150 }
151};
152} // namespace fims_distributions
153#endif
Definition fims_vector.hpp:27
size_type size() const
Returns the number of elements.
Definition fims_vector.hpp:273
Declares the DensityComponentBase class, which is the base class for all distribution functors.
#define FIMS_SIMULATE_F(F)
TMB macro that simulates data.
Definition interface.hpp:70
Base class for all module_name functors.
Definition density_components_base.hpp:32
fims::Vector< Type > lpdf_vec
Vector storing observation-level log-likelihood contributions.
Definition density_components_base.hpp:190
fims::Vector< Type > observed_values
Input value of distribution function for priors or random effects.
Definition density_components_base.hpp:61
bool simulate_flag
Boolean; if true, data are simulated from the distribution.
Definition density_components_base.hpp:200
Type & get_expected(size_t i)
Retrieve one expected value based on input_type and use_mean.
Definition density_components_base.hpp:122
Type lpdf
Total log probability density contribution of the distribution.
Definition density_components_base.hpp:180
size_t get_n_expected()
Get length of the active expected input vector.
Definition density_components_base.hpp:155
size_t get_n_x()
Get length of the active observed input vector.
Definition density_components_base.hpp:138
Type & get_observed(size_t i)
Retrieve one observed value based on input_type.
Definition density_components_base.hpp:83
fims::Vector< Type > * re
Pointer to random effects vector.
Definition density_components_base.hpp:47
std::string input_type
Classification of the input pathway for this distribution object. Options used by accessor methods ar...
Definition density_components_base.hpp:38
std::shared_ptr< fims_data_object::DataObject< Type > > data_observed_values
Observed data.
Definition density_components_base.hpp:41
std::vector< fims::Vector< Type > * > priors
Vector of pointers where each entry points to a prior parameter.
Definition density_components_base.hpp:56
Implements the LogNormalLPDF distribution functor used by FIMS to evaluate observation-level and tota...
Definition lognormal_lpdf.hpp:35
virtual ~LogNormalLPDF()
Destructor.
Definition lognormal_lpdf.hpp:50
LogNormalLPDF()
Constructor.
Definition lognormal_lpdf.hpp:46
fims::Vector< Type > log_sd
Natural log of the standard deviation of the distribution on the log scale. The argument can be a vec...
Definition lognormal_lpdf.hpp:42
virtual const Type evaluate()
Evaluates the lognormal log probability density function.
Definition lognormal_lpdf.hpp:63