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