FIMS  v0.8.1
Loading...
Searching...
No Matches
rcpp_distribution.hpp
Go to the documentation of this file.
1
9#ifndef FIMS_INTERFACE_RCPP_RCPP_OBJECTS_RCPP_DISTRIBUTION_HPP
10#define FIMS_INTERFACE_RCPP_RCPP_OBJECTS_RCPP_DISTRIBUTION_HPP
11
12#include "../../../distributions/distributions.hpp"
13#include "../../interface.hpp"
15
21 public:
25 static uint32_t id_g;
33 std::shared_ptr<std::vector<uint32_t>> key_m;
61 SharedString use_mean_m = fims::to_string("no");
67 static std::map<uint32_t, std::shared_ptr<DistributionsInterfaceBase>>
73
77 double lpdf_value = 0;
82 this->key_m = std::make_shared<std::vector<uint32_t>>();
84 /* Create instance of map: key is id and value is pointer to
85 DistributionsInterfaceBase */
86 // DistributionsInterfaceBase::live_objects[this->id_m] = this;
87 }
88
100
105
109 virtual uint32_t get_id() = 0;
110
119 virtual bool set_distribution_links(std::string input_type,
120 Rcpp::IntegerVector ids) {
121 return false;
122 }
123
147 virtual bool set_distribution_mean(double input_value) { return false; }
148
155 virtual bool set_observed_data(int observed_data_id) { return false; }
156
161 virtual double evaluate() = 0;
162};
163// static id of the DistributionsInterfaceBase object
165// local id of the DistributionsInterfaceBase object map relating the ID of the
166// DistributionsInterfaceBase to the DistributionsInterfaceBase objects
167std::map<uint32_t, std::shared_ptr<DistributionsInterfaceBase>>
169
175 public:
210
223
228
233 virtual uint32_t get_id() { return this->id_m; }
234
240 this->interface_observed_data_id_m.set(observed_data_id);
241 return true;
242 }
243
247 virtual bool set_distribution_mean(double input_value) {
248 this->expected_mean[0].initial_value_m = input_value;
249 this->expected_mean[0].estimation_type_m.set("fixed_effects");
250 this->use_mean_m.set(fims::to_string("yes"));
251 return true;
252 }
253
257 virtual bool set_distribution_links(std::string input_type,
258 Rcpp::IntegerVector ids) {
259 this->input_type_m.set(input_type);
260 this->key_m->resize(ids.size());
261 for (R_xlen_t i = 0; i < ids.size(); i++) {
262 this->key_m->at(i) = ids[i];
263 }
264 return true;
265 }
266
273 virtual double evaluate() {
275 dnorm.observed_values.resize(this->observed_values.size());
276 dnorm.expected_values.resize(this->expected_values.size());
277 dnorm.log_sd.resize(this->log_sd.size());
278 dnorm.expected_mean.resize(this->expected_mean.size());
279 for (size_t i = 0; i < this->observed_values.size(); i++) {
280 dnorm.observed_values[i] = this->observed_values[i].initial_value_m;
281 }
282 for (size_t i = 0; i < this->expected_values.size(); i++) {
283 dnorm.expected_values[i] = this->expected_values[i].initial_value_m;
284 }
285 for (size_t i = 0; i < this->log_sd.size(); i++) {
286 dnorm.log_sd[i] = this->log_sd[i].initial_value_m;
287 }
288 for (size_t i = 0; i < this->expected_mean.size(); i++) {
289 dnorm.expected_mean[i] = this->expected_mean[i].initial_value_m;
290 }
291 dnorm.use_mean = this->use_mean_m;
292 return dnorm.evaluate();
293 }
294
299 virtual void finalize() {
300 if (this->finalized) {
301 // log warning that finalize has been called more than once.
302 FIMS_WARNING_LOG("DnormDistribution " + fims::to_string(this->id_m) +
303 " has been finalized already.");
304 }
305
306 this->finalized = true; // indicate this has been called already
307
308 std::shared_ptr<fims_info::Information<double>> info =
310
312
313 // search for density component in Information
314 it = info->density_components.find(this->id_m);
315 // if not found, just return
316 if (it == info->density_components.end()) {
317 FIMS_WARNING_LOG("DnormDistribution " + fims::to_string(this->id_m) +
318 " not found in Information.");
319 return;
320 } else {
321 std::shared_ptr<fims_distributions::NormalLPDF<double>> dnorm =
322 std::dynamic_pointer_cast<fims_distributions::NormalLPDF<double>>(
323 it->second);
324
325 this->lpdf_value = dnorm->lpdf;
326
327 size_t n_x = dnorm->get_n_x();
328
329 // If input log_sd is a scalar, resize to n_x and fill with the scalar
330 // value
331 if (this->log_sd.size() != n_x) {
332 // If log_sd size == 1 (scalar), repeat the entry
333 if (this->log_sd.size() == 1) {
334 auto tmp = this->log_sd[0]; // copy the one log_sd param
335 this->log_sd.resize(n_x);
336 for (size_t i = 0; i < n_x; ++i) {
337 this->log_sd[i] = tmp; // copies all fields in Param
338 }
339 } else {
340 // Handle error
342 "log_sd size does not match number of observations and is not "
343 "scalar.");
344 }
345 }
346 for (size_t i = 0; i < n_x; i++) {
347 if (this->log_sd[i].estimation_type_m.get() == "constant") {
348 this->log_sd[i].final_value_m = this->log_sd[i].initial_value_m;
349 } else {
350 this->log_sd[i].final_value_m = dnorm->log_sd.get_force_scalar(i);
351 }
352 }
353
354 for (size_t i = 0; i < this->expected_mean.size(); i++) {
355 if (this->expected_mean[i].estimation_type_m.get() == "constant") {
356 this->expected_mean[i].final_value_m =
357 this->expected_mean[i].initial_value_m;
358 } else {
359 this->expected_mean[i].final_value_m = dnorm->expected_mean[i];
360 }
361 }
362
363 this->lpdf_vec = RealVector(n_x);
364 if (this->expected_values.size() == 1) {
365 this->expected_values.resize(n_x);
366 }
367 if (this->observed_values.size() == 1) {
368 this->observed_values.resize(n_x);
369 }
370
371 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
372 this->lpdf_vec[i] = dnorm->lpdf_vec[i];
373 this->expected_values[i].final_value_m = dnorm->get_expected(i);
374 this->observed_values[i].final_value_m = dnorm->get_observed(i);
375 }
376 }
377 }
378
386 virtual std::string to_json() {
387 std::stringstream ss;
388
389 ss << "{\n";
390 ss << " \"module_name\": \"density\",\n";
391 ss << " \"module_id\": " << this->id_m << ",\n";
392 ss << " \"module_type\": \"normal\",\n";
393 ss << " \"observed_data_id\" : " << this->interface_observed_data_id_m
394 << ",\n";
395 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
396 ss << " \"density_component\": {\n";
397 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
398 ss << " \"value\":[";
399 if (this->lpdf_vec.size() == 0) {
400 ss << "],\n";
401 } else {
402 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
403 ss << this->value_to_string(this->lpdf_vec[i]);
404 ss << ", ";
405 }
406 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
407
408 ss << "],\n";
409 }
410 ss << " \"expected_values\":[";
411 if (this->expected_values.size() == 0) {
412 ss << "],\n";
413 } else {
414 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
415 ss << this->value_to_string(this->expected_values[i].final_value_m)
416 << ", ";
417 }
418 ss << this->value_to_string(
419 this->expected_values[this->expected_values.size() - 1]
420 .final_value_m);
421 ss << "],\n";
422 }
423 ss << " \"log_sd_values\":[";
424 if (this->log_sd.size() == 0) {
425 ss << "],\n";
426 } else {
427 for (R_xlen_t i = 0; i < static_cast<R_xlen_t>(this->log_sd.size()) - 1;
428 i++) {
429 ss << this->value_to_string(this->log_sd[i].final_value_m) << ", ";
430 }
431 ss << this->value_to_string(
432 this->log_sd[this->log_sd.size() - 1].final_value_m)
433 << "],\n";
434 }
435 ss << " \"observed_values\":[";
436 if (this->observed_values.size() == 0) {
437 ss << "]\n";
438 } else {
439 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
440 ss << this->observed_values[i].final_value_m << ", ";
441 }
442 ss << this->observed_values[this->observed_values.size() - 1]
443 .final_value_m
444 << "]\n";
445 }
446 ss << " }}\n";
447 return ss.str();
448 }
449
450#ifdef TMB_MODEL
451
452 template <typename Type>
454 std::shared_ptr<fims_info::Information<Type>> info =
456
457 std::shared_ptr<fims_distributions::NormalLPDF<Type>> distribution =
458 std::make_shared<fims_distributions::NormalLPDF<Type>>();
459
460 // interface to data/parameter value
461
462 distribution->observed_data_id_m = interface_observed_data_id_m;
463 std::stringstream ss;
464 distribution->input_type = this->input_type_m;
465 distribution->key.resize(this->key_m->size());
466 for (size_t i = 0; i < this->key_m->size(); i++) {
467 distribution->key[i] = this->key_m->at(i);
468 }
469 distribution->id = this->id_m;
470 distribution->observed_values.resize(this->observed_values.size());
471 for (size_t i = 0; i < this->observed_values.size(); i++) {
472 distribution->observed_values[i] =
473 this->observed_values[i].initial_value_m;
474 }
475 // set relative info
476 distribution->expected_values.resize(this->expected_values.size());
477 for (size_t i = 0; i < this->expected_values.size(); i++) {
478 distribution->expected_values[i] =
479 this->expected_values[i].initial_value_m;
480 }
481 distribution->log_sd.resize(this->log_sd.size());
482 for (size_t i = 0; i < this->log_sd.size(); i++) {
483 distribution->log_sd[i] = this->log_sd[i].initial_value_m;
484 if (this->log_sd[i].estimation_type_m.get() == "fixed_effects") {
485 ss.str("");
486 ss << "dnorm." << this->id_m << ".log_sd." << this->log_sd[i].id_m;
487 info->RegisterParameterName(ss.str());
488 info->RegisterParameter(distribution->log_sd[i]);
489 }
490 if (this->log_sd[i].estimation_type_m.get() == "random_effects") {
491 FIMS_ERROR_LOG("standard deviations cannot be set to random effects");
492 }
493 }
494 info->variable_map[this->log_sd.id_m] = &(distribution)->log_sd;
495
496 distribution->use_mean = this->use_mean_m.get();
497 distribution->expected_mean.resize(this->expected_mean.size());
498 for (size_t i = 0; i < this->expected_mean.size(); i++) {
499 distribution->expected_mean[i] = this->expected_mean[i].initial_value_m;
500 if (this->expected_mean[i].estimation_type_m.get() == "fixed_effects") {
501 ss.str("");
502 ss << "dnorm." << this->id_m << ".expected_mean."
503 << this->expected_mean[i].id_m;
504 info->RegisterParameterName(ss.str());
505 info->RegisterParameter(distribution->expected_mean[i]);
506 }
507 if (this->expected_mean[i].estimation_type_m.get() == "random_effects") {
508 FIMS_ERROR_LOG("expected_mean cannot be set to random effects");
509 }
510 }
511 info->variable_map[this->expected_mean.id_m] =
513
514 info->density_components[distribution->id] = distribution;
515
516 return true;
517 }
518
523 virtual bool add_to_fims_tmb() {
526
527 return true;
528 }
529
530#endif
531};
532
538 public:
571
583
588
593 virtual uint32_t get_id() { return this->id_m; }
594
600 this->interface_observed_data_id_m.set(observed_data_id);
601 return true;
602 }
603
612 virtual bool set_distribution_links(std::string input_type,
613 Rcpp::IntegerVector ids) {
614 this->input_type_m.set(input_type);
615 this->key_m->resize(ids.size());
616 for (R_xlen_t i = 0; i < ids.size(); i++) {
617 this->key_m->at(i) = ids[i];
618 }
619 return true;
620 }
621
628 virtual double evaluate() {
630 dlnorm.observed_values.resize(this->observed_values.size());
631 dlnorm.expected_values.resize(this->expected_values.size());
632 dlnorm.log_sd.resize(this->log_sd.size());
633 // dlnorm.input_type = "prior";
634 for (size_t i = 0; i < this->observed_values.size(); i++) {
635 dlnorm.observed_values[i] = this->observed_values[i].initial_value_m;
636 }
637 for (size_t i = 0; i < this->expected_values.size(); i++) {
638 dlnorm.expected_values[i] = this->expected_values[i].initial_value_m;
639 }
640 for (size_t i = 0; i < this->log_sd.size(); i++) {
641 dlnorm.log_sd[i] = this->log_sd[i].initial_value_m;
642 }
643 return dlnorm.evaluate();
644 }
645
650 virtual void finalize() {
651 if (this->finalized) {
652 // log warning that finalize has been called more than once.
653 FIMS_WARNING_LOG("LogNormalLPDF " + fims::to_string(this->id_m) +
654 " has been finalized already.");
655 }
656
657 this->finalized = true; // indicate this has been called already
658
659 std::shared_ptr<fims_info::Information<double>> info =
661
663
664 // search for density component in Information
665 it = info->density_components.find(this->id_m);
666 // if not found, just return
667 if (it == info->density_components.end()) {
668 FIMS_WARNING_LOG("LogNormalLPDF " + fims::to_string(this->id_m) +
669 " not found in Information.");
670 return;
671 } else {
672 std::shared_ptr<fims_distributions::LogNormalLPDF<double>> dlnorm =
673 std::dynamic_pointer_cast<fims_distributions::LogNormalLPDF<double>>(
674 it->second);
675
676 this->lpdf_value = dlnorm->lpdf;
677
678 size_t n_x = dlnorm->get_n_x();
679
680 if (this->log_sd.size() != n_x) {
681 // If log_sd size == 1 (scalar), repeat the entry
682 if (this->log_sd.size() == 1) {
683 auto tmp = this->log_sd[0]; // copy the one log_sd param
684 this->log_sd.resize(n_x);
685 for (size_t i = 0; i < n_x; ++i) {
686 this->log_sd[i] = tmp; // copies all fields in Param
687 }
688 } else {
689 // Handle error
691 "log_sd size does not match number of observations and is not "
692 "scalar.");
693 }
694 }
695
696 for (size_t i = 0; i < n_x; i++) {
697 if (this->log_sd[i].estimation_type_m.get() == "constant") {
698 this->log_sd[i].final_value_m = this->log_sd[i].initial_value_m;
699 } else {
700 this->log_sd[i].final_value_m = dlnorm->log_sd.get_force_scalar(i);
701 }
702 }
703
704 this->lpdf_vec = RealVector(n_x);
705 if (this->expected_values.size() == 1) {
706 this->expected_values.resize(n_x);
707 }
708 if (this->observed_values.size() == 1) {
709 this->observed_values.resize(n_x);
710 }
711 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
712 this->lpdf_vec[i] = dlnorm->lpdf_vec[i];
713 this->expected_values[i].final_value_m = dlnorm->get_expected(i);
714 this->observed_values[i].final_value_m = dlnorm->get_observed(i);
715 }
716 }
717 }
718
726 virtual std::string to_json() {
727 std::stringstream ss;
728
729 ss << "{\n";
730 ss << " \"module_name\": \"density\",\n";
731 ss << " \"module_id\": " << this->id_m << ",\n";
732 ss << " \"module_type\": \"log_normal\",\n";
733 ss << " \"observed_data_id\" : " << this->interface_observed_data_id_m
734 << ",\n";
735 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
736 ss << " \"density_component\": {\n";
737 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
738 ss << " \"value\":[";
739 if (this->lpdf_vec.size() == 0) {
740 ss << "],\n";
741 } else {
742 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
743 ss << this->value_to_string(this->lpdf_vec[i]) << ", ";
744 }
745 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
746
747 ss << "],\n";
748 }
749 ss << " \"expected_values\":[";
750 if (this->expected_values.size() == 0) {
751 ss << "],\n";
752 } else {
753 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
754 ss << this->value_to_string(this->expected_values[i].final_value_m)
755 << ", ";
756 }
757 ss << this->value_to_string(
758 this->expected_values[this->expected_values.size() - 1]
759 .final_value_m);
760
761 ss << "],\n";
762 }
763 ss << " \"log_sd_values\":[";
764 if (this->log_sd.size() == 0) {
765 ss << "],\n";
766 } else {
767 for (R_xlen_t i = 0; i < static_cast<R_xlen_t>(this->log_sd.size()) - 1;
768 i++) {
769 ss << this->value_to_string(this->log_sd[i].final_value_m) << ", ";
770 }
771 ss << this->value_to_string(
772 this->log_sd[this->log_sd.size() - 1].final_value_m)
773 << "],\n";
774 }
775 ss << " \"observed_values\":[";
776 if (this->observed_values.size() == 0) {
777 ss << "]\n";
778 } else {
779 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
780 ss << this->observed_values[i].final_value_m << ", ";
781 }
782 ss << this->observed_values[this->observed_values.size() - 1]
783 .final_value_m
784 << "]\n";
785 }
786 ss << " }}\n";
787 return ss.str();
788 }
789
790#ifdef TMB_MODEL
791
792 template <typename Type>
794 std::shared_ptr<fims_info::Information<Type>> info =
796
797 std::shared_ptr<fims_distributions::LogNormalLPDF<Type>> distribution =
798 std::make_shared<fims_distributions::LogNormalLPDF<Type>>();
799
800 // set relative info
801 distribution->id = this->id_m;
802 std::stringstream ss;
803 distribution->observed_data_id_m = interface_observed_data_id_m;
804 distribution->input_type = this->input_type_m;
805 distribution->key.resize(this->key_m->size());
806 for (size_t i = 0; i < this->key_m->size(); i++) {
807 distribution->key[i] = this->key_m->at(i);
808 }
809 distribution->observed_values.resize(this->observed_values.size());
810 for (size_t i = 0; i < this->observed_values.size(); i++) {
811 distribution->observed_values[i] =
812 this->observed_values[i].initial_value_m;
813 }
814 // set relative info
815 distribution->expected_values.resize(this->expected_values.size());
816 for (size_t i = 0; i < this->expected_values.size(); i++) {
817 distribution->expected_values[i] =
818 this->expected_values[i].initial_value_m;
819 }
820 distribution->log_sd.resize(this->log_sd.size());
821 for (size_t i = 0; i < this->log_sd.size(); i++) {
822 distribution->log_sd[i] = this->log_sd[i].initial_value_m;
823 if (this->log_sd[i].estimation_type_m.get() == "fixed_effects") {
824 ss.str("");
825 ss << "dlnorm." << this->id_m << ".log_sd." << this->log_sd[i].id_m;
826 info->RegisterParameterName(ss.str());
827 info->RegisterParameter(distribution->log_sd[i]);
828 }
829 if (this->log_sd[i].estimation_type_m.get() == "random_effects") {
830 FIMS_ERROR_LOG("standard deviations cannot be set to random effects");
831 }
832 }
833 info->variable_map[this->log_sd.id_m] = &(distribution)->log_sd;
834
835 info->density_components[distribution->id] = distribution;
836
837 return true;
838 }
839
844 virtual bool add_to_fims_tmb() {
847
848 return true;
849 }
850
851#endif
852};
853
859 public:
886
896
909
918 virtual uint32_t get_id() { return this->id_m; }
919
925 this->interface_observed_data_id_m.set(observed_data_id);
926 return true;
927 }
928
937 virtual bool set_distribution_links(std::string input_type,
938 Rcpp::IntegerVector ids) {
939 this->input_type_m.set(input_type);
940 this->key_m->resize(ids.size());
941 for (R_xlen_t i = 0; i < ids.size(); i++) {
942 this->key_m->at(i) = ids[i];
943 }
944 return true;
945 }
946
952 void set_note(std::string note) { this->notes.set(note); }
953
959 virtual double evaluate() {
961 // Declare TMBVector in this scope
962 dmultinom.observed_values.resize(this->observed_values.size());
963 dmultinom.expected_values.resize(this->expected_values.size());
964 for (size_t i = 0; i < observed_values.size(); i++) {
965 dmultinom.observed_values[i] = this->observed_values[i].initial_value_m;
966 }
967 for (size_t i = 0; i < expected_values.size(); i++) {
968 dmultinom.expected_values[i] = this->expected_values[i].initial_value_m;
969 }
970 dmultinom.dims.resize(2);
971 dmultinom.dims[0] = this->dims[0];
972 dmultinom.dims[1] = this->dims[1];
973 return dmultinom.evaluate();
974 }
975
976 void finalize() {
977 if (this->finalized) {
978 // log warning that finalize has been called more than once.
979 FIMS_WARNING_LOG("DmultinomDistributions " +
980 fims::to_string(this->id_m) +
981 " has been finalized already.");
982 }
983
984 this->finalized = true; // indicate this has been called already
985
986 std::shared_ptr<fims_info::Information<double>> info =
988
990
991 // search for density component in Information
992 it = info->density_components.find(this->id_m);
993 // if not found, just return
994 if (it == info->density_components.end()) {
995 FIMS_WARNING_LOG("DmultinomDistributions " + fims::to_string(this->id_m) +
996 " not found in Information.");
997 return;
998 } else {
999 std::shared_ptr<fims_distributions::MultinomialLPMF<double>> dmultinom =
1000 std::dynamic_pointer_cast<
1002
1003 this->lpdf_value = dmultinom->lpdf;
1004
1005 size_t n_x = dmultinom->lpdf_vec.size();
1006 this->lpdf_vec = Rcpp::NumericVector(n_x);
1007 if (this->expected_values.size() != n_x) {
1008 this->expected_values.resize(n_x);
1009 }
1010 if (this->observed_values.size() != n_x) {
1011 this->observed_values.resize(n_x);
1012 }
1013 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
1014 this->lpdf_vec[i] = dmultinom->lpdf_vec[i];
1015 this->expected_values[i].final_value_m = dmultinom->get_expected(i);
1016 if (dmultinom->input_type != "data") {
1017 this->observed_values[i].final_value_m = dmultinom->get_observed(i);
1018 }
1019 }
1020 if (dmultinom->input_type == "data") {
1021 dims.resize(2);
1022 dims[0] = dmultinom->dims[0];
1023 dims[1] = dmultinom->dims[1];
1024 for (size_t i = 0; i < dims[0]; i++) {
1025 for (size_t j = 0; j < dims[1]; j++) {
1026 size_t idx = (i * dims[1]) + j;
1027 this->observed_values[idx].final_value_m = dmultinom->get_observed(
1028 static_cast<size_t>(i), static_cast<size_t>(j));
1029 }
1030 }
1031 }
1032 }
1033 }
1034
1042 virtual std::string to_json() {
1043 std::stringstream ss;
1044
1045 ss << "{\n";
1046 ss << " \"module_name\": \"density\",\n";
1047 ss << " \"module_id\": " << this->id_m << ",\n";
1048 ss << " \"module_type\": \"multinomial\",\n";
1049 ss << "\"observed_data_id\" : " << this->interface_observed_data_id_m
1050 << ",\n";
1051 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
1052 ss << " \"density_component\": {\n";
1053 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
1054 ss << " \"value\":[";
1055 if (this->lpdf_vec.size() == 0) {
1056 ss << "],\n";
1057 } else {
1058 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
1059 ss << this->value_to_string(this->lpdf_vec[i]);
1060 ss << ", ";
1061 }
1062 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
1063
1064 ss << "],\n";
1065 }
1066 ss << " \"expected_values\":[";
1067 if (this->expected_values.size() == 0) {
1068 ss << "],\n";
1069 } else {
1070 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
1071 ss << this->value_to_string(this->expected_values[i].final_value_m)
1072 << ", ";
1073 }
1074 ss << this->value_to_string(
1075 this->expected_values[this->expected_values.size() - 1]
1076 .final_value_m);
1077
1078 ss << "],\n";
1079 }
1080 // no log_sd_values for multinomial
1081 ss << " \"observed_values\":[";
1082 if (this->observed_values.size() == 0) {
1083 ss << "]\n";
1084 } else {
1085 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
1086 ss << this->observed_values[i].final_value_m << ", ";
1087 }
1088 ss << this->observed_values[this->observed_values.size() - 1]
1089 .final_value_m
1090 << "]\n";
1091 }
1092 ss << " }}\n";
1093 return ss.str();
1094 }
1095
1096#ifdef TMB_MODEL
1097
1098 template <typename Type>
1100 FIMS_INFO_LOG("Adding multinomial to FIMS.");
1101 std::shared_ptr<fims_info::Information<Type>> info =
1103
1104 std::shared_ptr<fims_distributions::MultinomialLPMF<Type>> distribution =
1105 std::make_shared<fims_distributions::MultinomialLPMF<Type>>();
1106
1107 distribution->id = this->id_m;
1108 distribution->observed_data_id_m = interface_observed_data_id_m;
1109 distribution->input_type = this->input_type_m;
1110 distribution->key.resize(this->key_m->size());
1111 for (size_t i = 0; i < this->key_m->size(); i++) {
1112 distribution->key[i] = this->key_m->at(i);
1113 }
1114 distribution->observed_values.resize(this->observed_values.size());
1115 for (size_t i = 0; i < this->observed_values.size(); i++) {
1116 distribution->observed_values[i] =
1117 this->observed_values[i].initial_value_m;
1118 }
1119 // set relative info
1120 distribution->expected_values.resize(this->expected_values.size());
1121 for (size_t i = 0; i < this->expected_values.size(); i++) {
1122 distribution->expected_values[i] =
1123 this->expected_values[i].initial_value_m;
1124 }
1125
1126 info->density_components[distribution->id] = distribution;
1127 FIMS_INFO_LOG("Done adding multinomial to FIMS.");
1128 return true;
1129 }
1130
1131 virtual bool add_to_fims_tmb() {
1134
1135 return true;
1136 }
1137
1138#endif
1139};
1140
1141#endif
Rcpp interface that serves as the parent class for Rcpp distribution interfaces. This type should be ...
Definition rcpp_distribution.hpp:20
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:155
double lpdf_value
The log probability density function value.
Definition rcpp_distribution.hpp:77
virtual ~DistributionsInterfaceBase()
The destructor.
Definition rcpp_distribution.hpp:104
virtual double evaluate()=0
A method for each child distribution interface object to inherit so each distribution can have an eva...
std::shared_ptr< std::vector< uint32_t > > key_m
The unique ID for the variable map that points to a fims::Vector.
Definition rcpp_distribution.hpp:33
DistributionsInterfaceBase(const DistributionsInterfaceBase &other)
Construct a new Distributions Interface Base object.
Definition rcpp_distribution.hpp:94
uint32_t id_m
The local ID of the DistributionsInterfaceBase object.
Definition rcpp_distribution.hpp:29
virtual uint32_t get_id()=0
Get the ID for the child distribution interface objects to inherit.
SharedString use_mean_m
Control flag indicating whether to use the expected mean in the distribution calculations.
Definition rcpp_distribution.hpp:61
SharedInt interface_observed_data_id_m
The ID of the observed data object, which is set to -999.
Definition rcpp_distribution.hpp:72
virtual bool set_distribution_mean(double input_value)
Set the expected mean value for the distribution.
Definition rcpp_distribution.hpp:147
virtual bool set_distribution_links(std::string input_type, Rcpp::IntegerVector ids)
Sets pointers for data observations, random effects, or priors.
Definition rcpp_distribution.hpp:119
static uint32_t id_g
The static ID of the DistributionsInterfaceBase object.
Definition rcpp_distribution.hpp:25
static std::map< uint32_t, std::shared_ptr< DistributionsInterfaceBase > > live_objects
The map associating the ID of the DistributionsInterfaceBase to the DistributionsInterfaceBase object...
Definition rcpp_distribution.hpp:68
SharedString input_type_m
The type of density input. The options are prior, re, or data.
Definition rcpp_distribution.hpp:37
DistributionsInterfaceBase()
The constructor.
Definition rcpp_distribution.hpp:81
The Rcpp interface for Dlnorm to instantiate from R: dlnorm_ <- methods::new(DlnormDistribution).
Definition rcpp_distribution.hpp:537
virtual bool set_distribution_links(std::string input_type, Rcpp::IntegerVector ids)
Sets pointers for data observations, random effects, or priors.
Definition rcpp_distribution.hpp:612
DlnormDistributionsInterface(const DlnormDistributionsInterface &other)
Construct a new Dlnorm Distributions Interface object.
Definition rcpp_distribution.hpp:577
DlnormDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:565
ParameterVector expected_values
The expected values, which would be the mean of log(x) for this distribution.
Definition rcpp_distribution.hpp:547
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:726
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:593
virtual double evaluate()
Evaluate lognormal probability density function (pdf). The natural log of the pdf is returned.
Definition rcpp_distribution.hpp:628
ParameterVector observed_values
Observed data.
Definition rcpp_distribution.hpp:542
ParameterVector log_sd
The uncertainty, which would be the natural logarithm of the standard deviation (sd) of log(x) for th...
Definition rcpp_distribution.hpp:555
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:599
virtual void finalize()
Extracts the derived quantities from Information to the Rcpp object.
Definition rcpp_distribution.hpp:650
virtual ~DlnormDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:587
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:560
The Rcpp interface for Dmultinom to instantiate from R: dmultinom_ <- methods::new(DmultinomDistribut...
Definition rcpp_distribution.hpp:858
virtual bool set_distribution_links(std::string input_type, Rcpp::IntegerVector ids)
Sets pointers for data observations, random effects, or priors.
Definition rcpp_distribution.hpp:937
SharedString notes
TODO: document this.
Definition rcpp_distribution.hpp:885
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:918
ParameterVector expected_values
The expected values, which should be a vector of length K where each value specifies the probability ...
Definition rcpp_distribution.hpp:869
void finalize()
Extracts derived quantities back to the Rcpp interface object from the Information object.
Definition rcpp_distribution.hpp:976
virtual double evaluate()
Definition rcpp_distribution.hpp:959
ParameterVector observed_values
Observed data, which should be a vector of length K of integers.
Definition rcpp_distribution.hpp:863
void set_note(std::string note)
Set the note object.
Definition rcpp_distribution.hpp:952
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:924
RealVector dims
The dimensions of the number of rows and columns of the multivariate dataset.
Definition rcpp_distribution.hpp:874
DmultinomDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:890
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:1042
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:879
DmultinomDistributionsInterface(const DmultinomDistributionsInterface &other)
Construct a new Dmultinom Distributions Interface object.
Definition rcpp_distribution.hpp:902
virtual ~DmultinomDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:913
The Rcpp interface for Dnorm to instantiate from R: dnorm_ <- methods::new(DnormDistribution).
Definition rcpp_distribution.hpp:174
DnormDistributionsInterface(const DnormDistributionsInterface &other)
Construct a new Dnorm Distributions Interface object.
Definition rcpp_distribution.hpp:216
ParameterVector expected_mean
The expected mean, which would be the mean of x for this distribution.
Definition rcpp_distribution.hpp:189
DnormDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:204
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:386
virtual ~DnormDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:227
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:199
virtual double evaluate()
Evaluate normal probability density function (pdf). The natural log of the pdf is returned.
Definition rcpp_distribution.hpp:273
ParameterVector expected_values
The expected values, which would be the mean of x for this distribution.
Definition rcpp_distribution.hpp:184
virtual void finalize()
Extracts the derived quantities from Information to the Rcpp object.
Definition rcpp_distribution.hpp:299
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:239
ParameterVector observed_values
Observed data.
Definition rcpp_distribution.hpp:179
virtual bool set_distribution_mean(double input_value)
Set the expected mean value for the distribution.
Definition rcpp_distribution.hpp:247
ParameterVector log_sd
The uncertainty, which would be the standard deviation of x for the normal distribution.
Definition rcpp_distribution.hpp:194
virtual bool set_distribution_links(std::string input_type, Rcpp::IntegerVector ids)
Sets pointers for data observations, random effects, or priors.
Definition rcpp_distribution.hpp:257
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:233
Base class for all interface objects.
Definition rcpp_interface_base.hpp:574
bool finalized
Is the object already finalized? The default is false.
Definition rcpp_interface_base.hpp:579
std::string value_to_string(double value)
Report the parameter value as a string.
Definition rcpp_interface_base.hpp:615
static std::vector< std::shared_ptr< FIMSRcppInterfaceBase > > fims_interface_objects
FIMS interface object vectors.
Definition rcpp_interface_base.hpp:584
virtual bool add_to_fims_tmb()
A virtual method to inherit to add objects to the TMB model.
Definition rcpp_interface_base.hpp:589
An Rcpp interface class that defines the ParameterVector class.
Definition rcpp_interface_base.hpp:142
void resize(size_t size)
Resizes a ParameterVector to the desired length.
Definition rcpp_interface_base.hpp:287
size_t size()
Returns the size of a ParameterVector.
Definition rcpp_interface_base.hpp:280
uint32_t id_m
The local ID of the Parameter object.
Definition rcpp_interface_base.hpp:155
Parameter & get(size_t pos)
An internal accessor for calling a position of a ParameterVector from R.
Definition rcpp_interface_base.hpp:259
void set(size_t pos, const Parameter &p)
An internal setter for setting a position of a ParameterVector from R.
Definition rcpp_interface_base.hpp:275
An Rcpp interface class that defines the RealVector class.
Definition rcpp_interface_base.hpp:374
size_t size()
Returns the size of a RealVector.
Definition rcpp_interface_base.hpp:535
void resize(size_t size)
Resizes a RealVector to the desired length.
Definition rcpp_interface_base.hpp:542
A class that provides shared ownership of an integer value.
Definition rcpp_shared_primitive.hpp:24
void set(int val)
Change the value of the integer.
Definition rcpp_shared_primitive.hpp:129
A class that provides shared ownership of a string.
Definition rcpp_shared_primitive.hpp:1508
std::string get() const
Retrieves the string value managed by the object.
Definition rcpp_shared_primitive.hpp:1612
void set(std::string val)
Sets the string value of the object.
Definition rcpp_shared_primitive.hpp:1621
std::map< uint32_t, std::shared_ptr< fims_distributions::DensityComponentBase< Type > > > density_components
Definition information.hpp:126
std::map< uint32_t, std::shared_ptr< fims_distributions::DensityComponentBase< Type > > >::iterator density_components_iterator
Definition information.hpp:131
static std::shared_ptr< Information< Type > > GetInstance()
Returns a singleton Information object for type T.
Definition information.hpp:238
#define FIMS_WARNING_LOG(MESSAGE)
Definition def.hpp:545
#define FIMS_INFO_LOG(MESSAGE)
Definition def.hpp:537
#define FIMS_ERROR_LOG(MESSAGE)
Definition def.hpp:553
void clear_internal()
Clears the internal objects.
Definition rcpp_interface.hpp:239
The Rcpp interface to declare objects that are used ubiquitously throughout the Rcpp interface,...
fims::Vector< Type > observed_values
Input value of distribution function for priors or random effects.
Definition density_components_base.hpp:61
Type lpdf
Total log probability density contribution of the distribution.
Definition density_components_base.hpp:180
Implements the LogNormalLPDF distribution functor used by FIMS to evaluate observation-level and tota...
Definition lognormal_lpdf.hpp:35
Implements the MultinomialLPMF distribution functor used by FIMS to evaluate the observation-level an...
Definition multinomial_lpmf.hpp:37
Implements the NormalLPDF distribution functor used by FIMS to evaluate observation-level and total l...
Definition normal_lpdf.hpp:33
virtual const Type evaluate()
Evaluates the normal log probability density function.
Definition normal_lpdf.hpp:60