FIMS  v0.9.2
Loading...
Searching...
No Matches
normal_lpdf.hpp
Go to the documentation of this file.
1
11#ifndef NORMAL_LPDF
12#define NORMAL_LPDF
13
14#include "../../common/def.hpp"
16#include "../../common/fims_vector.hpp"
17namespace fims_distributions {
31template <typename Type>
32struct NormalLPDF : public DensityComponentBase<Type> {
40
44
47 virtual ~NormalLPDF() {}
48
59 virtual const Type evaluate() {
60 // set vector size based on input type (prior, process, or data)
61 size_t n_x = this->get_n_x();
62 // get expected value vector size
63 size_t n_expected = this->get_n_expected();
64 // setup vector for recording the log probability density function values
65 this->lpdf_vec.resize(n_x);
66 std::fill(this->lpdf_vec.begin(), this->lpdf_vec.end(),
67 static_cast<Type>(0));
68 this->lpdf = static_cast<Type>(0);
69
70 // Dimension checks
71 if (n_x != n_expected) {
72 if (n_expected == 1) {
73 n_expected = n_x;
74 } else if (n_x > n_expected) {
75 n_x = n_expected;
76 }
77 }
78
79 if (this->log_sd.size() > 1 && n_x != this->log_sd.size()) {
80 throw std::invalid_argument(
81 "NormalLPDF::Vector index out of bounds. The size of observed data "
82 "does not equal the size of the log_sd vector. The observed data "
83 "vector is of size " +
84 fims::to_string(n_x) + " and the log_sd vector is of size " +
85 fims::to_string(this->log_sd.size()));
86 }
87
88 for (size_t i = 0; i < n_x; i++) {
89#ifdef TMB_MODEL
90 if (this->input_type == "data") {
91 // if data, check if there are any NA values and skip lpdf calculation
92 // if there are
93 if (this->get_observed(i) != this->data_observed_values->na_value) {
94 this->lpdf_vec[i] =
95 dnorm(this->get_observed(i), this->get_expected(i),
96 fims_math::exp(log_sd.get_force_scalar(i)), true);
97 } else {
98 this->lpdf_vec[i] = 0;
99 }
100 // if not data (i.e. prior or process), use x vector instead of
101 // observed_values
102 } else {
103 this->lpdf_vec[i] =
104 dnorm(this->get_observed(i), this->get_expected(i),
105 fims_math::exp(log_sd.get_force_scalar(i)), true);
106 }
107 this->lpdf += this->lpdf_vec[i];
108 if (this->simulate_flag) {
109 FIMS_SIMULATE_F(this->of) {
110 if (this->input_type == "data") {
111 this->data_observed_values->at(i) =
112 rnorm(this->get_expected(i),
113 fims_math::exp(log_sd.get_force_scalar(i)));
114 }
115 if (this->input_type == "random_effects") {
116 (*this->re)[i] = rnorm(this->get_expected(i),
117 fims_math::exp(log_sd.get_force_scalar(i)));
118 }
119 if (this->input_type == "prior") {
120 (*(this->priors[i]))[0] =
121 rnorm(this->get_expected(i),
122 fims_math::exp(log_sd.get_force_scalar(i)));
123 }
124 }
125 }
126#endif
127 /* osa not working yet
128 if(osa_flag){//data observation type implements osa residuals
129 //code for osa cdf method
130 this->lpdf_vec[i] = this->keep.cdf_lower[i] * log(
131 pnorm(this->observed_values[i], this->get_expected(i), sd[i]) );
132 this->lpdf_vec[i] = this->keep.cdf_upper[i] * log( 1.0 -
133 pnorm(this->observed_values[i], this->get_expected(i), sd[i]) );
134 } */
135 }
136#ifdef TMB_MODEL
137 vector<Type> normal_observed_values = this->observed_values.to_tmb();
138#endif
139 return (this->lpdf);
140 }
141};
142
143} // namespace fims_distributions
144#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:71
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 NormalLPDF distribution functor used by FIMS to evaluate observation-level and total l...
Definition normal_lpdf.hpp:32
NormalLPDF()
Constructor.
Definition normal_lpdf.hpp:43
fims::Vector< Type > log_sd
The natural log of the standard deviation of the distribution. The argument can be a vector or scalar...
Definition normal_lpdf.hpp:39
virtual ~NormalLPDF()
Destructor.
Definition normal_lpdf.hpp:47
virtual const Type evaluate()
Evaluates the normal log probability density function.
Definition normal_lpdf.hpp:59