FIMS  v0.9.2
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
13#include "../../../distributions/distributions.hpp"
14
20 public:
24 static uint32_t id_g;
32 std::shared_ptr<std::vector<uint32_t>> key_m;
60 SharedString use_mean_m = fims::to_string("no");
66 static std::map<uint32_t, std::shared_ptr<DistributionsInterfaceBase>>
72
76 double lpdf_value = 0;
81 this->key_m = std::make_shared<std::vector<uint32_t>>();
83 /* Create instance of map: key is id and value is pointer to
84 DistributionsInterfaceBase */
85 // DistributionsInterfaceBase::live_objects[this->id_m] = this;
86 }
87
99
104
108 virtual uint32_t get_id() = 0;
109
118 virtual bool set_distribution_links(std::string input_type,
119 Rcpp::IntegerVector ids) {
120 return false;
121 }
122
146 virtual bool set_distribution_mean(double input_value) { return false; }
147
154 virtual bool set_observed_data(int observed_data_id) { return false; }
155
160 virtual double evaluate() = 0;
161};
162// static id of the DistributionsInterfaceBase object
164// local id of the DistributionsInterfaceBase object map relating the ID of the
165// DistributionsInterfaceBase to the DistributionsInterfaceBase objects
166std::map<uint32_t, std::shared_ptr<DistributionsInterfaceBase>>
168
174 public:
209
222
227
232 virtual uint32_t get_id() { return this->id_m; }
233
239 this->interface_observed_data_id_m.set(observed_data_id);
240 return true;
241 }
242
246 virtual bool set_distribution_mean(double input_value) {
247 this->expected_mean[0].initial_value_m = input_value;
248 this->expected_mean[0].estimation_type_m.set("fixed_effects");
249 this->use_mean_m.set(fims::to_string("yes"));
250 return true;
251 }
252
256 virtual bool set_distribution_links(std::string input_type,
257 Rcpp::IntegerVector ids) {
258 this->input_type_m.set(input_type);
259 this->key_m->resize(ids.size());
260 for (R_xlen_t i = 0; i < ids.size(); i++) {
261 this->key_m->at(i) = ids[i];
262 }
263 return true;
264 }
265
272 virtual double evaluate() {
274 dnorm.observed_values.resize(this->observed_values.size());
275 dnorm.expected_values.resize(this->expected_values.size());
276 dnorm.log_sd.resize(this->log_sd.size());
277 dnorm.expected_mean.resize(this->expected_mean.size());
278 for (size_t i = 0; i < this->observed_values.size(); i++) {
279 dnorm.observed_values[i] = this->observed_values[i].initial_value_m;
280 }
281 for (size_t i = 0; i < this->expected_values.size(); i++) {
282 dnorm.expected_values[i] = this->expected_values[i].initial_value_m;
283 }
284 for (size_t i = 0; i < this->log_sd.size(); i++) {
285 dnorm.log_sd[i] = this->log_sd[i].initial_value_m;
286 }
287 for (size_t i = 0; i < this->expected_mean.size(); i++) {
288 dnorm.expected_mean[i] = this->expected_mean[i].initial_value_m;
289 }
290 dnorm.use_mean = this->use_mean_m;
291 return dnorm.evaluate();
292 }
293
298 virtual void finalize() {
299 if (this->finalized) {
300 // log warning that finalize has been called more than once.
301 FIMS_WARNING_LOG("DnormDistribution " + fims::to_string(this->id_m) +
302 " has been finalized already.");
303 }
304
305 this->finalized = true; // indicate this has been called already
306
307 std::shared_ptr<fims_info::Information<double>> info =
309
311
312 // search for density component in Information
313 it = info->density_components.find(this->id_m);
314 // if not found, just return
315 if (it == info->density_components.end()) {
316 FIMS_WARNING_LOG("DnormDistribution " + fims::to_string(this->id_m) +
317 " not found in Information.");
318 return;
319 } else {
320 std::shared_ptr<fims_distributions::NormalLPDF<double>> dnorm =
321 std::dynamic_pointer_cast<fims_distributions::NormalLPDF<double>>(
322 it->second);
323
324 this->lpdf_value = dnorm->lpdf;
325
326 size_t n_x = dnorm->get_n_x();
327
328 // If input log_sd is a scalar, resize to n_x and fill with the scalar
329 // value
330 if (this->log_sd.size() != n_x) {
331 // If log_sd size == 1 (scalar), repeat the entry
332 if (this->log_sd.size() == 1) {
333 auto tmp = this->log_sd[0]; // copy the one log_sd param
334 this->log_sd.resize(n_x);
335 for (size_t i = 0; i < n_x; ++i) {
336 this->log_sd[i] = tmp; // copies all fields in Param
337 }
338 } else {
339 // Handle error
341 "log_sd size does not match number of observations and is not "
342 "scalar.");
343 }
344 }
345 for (size_t i = 0; i < n_x; i++) {
346 if (this->log_sd[i].estimation_type_m.get() == "constant") {
347 this->log_sd[i].final_value_m = this->log_sd[i].initial_value_m;
348 } else {
349 this->log_sd[i].final_value_m = dnorm->log_sd.get_force_scalar(i);
350 }
351 }
352
353 for (size_t i = 0; i < this->expected_mean.size(); i++) {
354 if (this->expected_mean[i].estimation_type_m.get() == "constant") {
355 this->expected_mean[i].final_value_m =
356 this->expected_mean[i].initial_value_m;
357 } else {
358 this->expected_mean[i].final_value_m = dnorm->expected_mean[i];
359 }
360 }
361
362 this->lpdf_vec = RealVector(n_x);
363 if (this->expected_values.size() == 1) {
364 this->expected_values.resize(n_x);
365 }
366 if (this->observed_values.size() == 1) {
367 this->observed_values.resize(n_x);
368 }
369
370 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
371 this->lpdf_vec[i] = dnorm->lpdf_vec[i];
372 this->expected_values[i].final_value_m = dnorm->get_expected(i);
373 this->observed_values[i].final_value_m = dnorm->get_observed(i);
374 }
375 }
376 }
377
385 virtual std::string to_json() {
386 std::stringstream ss;
387
388 ss << "{\n";
389 ss << " \"module_name\": \"density\",\n";
390 ss << " \"module_id\": " << this->id_m << ",\n";
391 ss << " \"module_type\": \"normal\",\n";
392 ss << " \"observed_data_id\" : " << this->interface_observed_data_id_m
393 << ",\n";
394 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
395 ss << " \"density_component\": {\n";
396 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
397 ss << " \"value\":[";
398 if (this->lpdf_vec.size() == 0) {
399 ss << "],\n";
400 } else {
401 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
402 ss << this->value_to_string(this->lpdf_vec[i]);
403 ss << ", ";
404 }
405 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
406
407 ss << "],\n";
408 }
409 ss << " \"expected_values\":[";
410 if (this->expected_values.size() == 0) {
411 ss << "],\n";
412 } else {
413 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
414 ss << this->value_to_string(this->expected_values[i].final_value_m)
415 << ", ";
416 }
417 ss << this->value_to_string(
418 this->expected_values[this->expected_values.size() - 1]
419 .final_value_m);
420 ss << "],\n";
421 }
422 ss << " \"log_sd_values\":[";
423 if (this->log_sd.size() == 0) {
424 ss << "],\n";
425 } else {
426 for (R_xlen_t i = 0; i < static_cast<R_xlen_t>(this->log_sd.size()) - 1;
427 i++) {
428 ss << this->value_to_string(this->log_sd[i].final_value_m) << ", ";
429 }
430 ss << this->value_to_string(
431 this->log_sd[this->log_sd.size() - 1].final_value_m)
432 << "],\n";
433 }
434 ss << " \"observed_values\":[";
435 if (this->observed_values.size() == 0) {
436 ss << "]\n";
437 } else {
438 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
439 ss << this->observed_values[i].final_value_m << ", ";
440 }
441 ss << this->observed_values[this->observed_values.size() - 1]
442 .final_value_m
443 << "]\n";
444 }
445 ss << " }}\n";
446 return ss.str();
447 }
448
449#ifdef TMB_MODEL
450
451 template <typename Type>
453 std::shared_ptr<fims_info::Information<Type>> info =
455
456 std::shared_ptr<fims_distributions::NormalLPDF<Type>> distribution =
457 std::make_shared<fims_distributions::NormalLPDF<Type>>();
458
459 // interface to data/parameter value
460
461 distribution->observed_data_id_m = interface_observed_data_id_m;
462 std::stringstream ss;
463 distribution->input_type = this->input_type_m;
464 distribution->key.resize(this->key_m->size());
465 for (size_t i = 0; i < this->key_m->size(); i++) {
466 distribution->key[i] = this->key_m->at(i);
467 }
468 distribution->id = this->id_m;
469 distribution->observed_values.resize(this->observed_values.size());
470 for (size_t i = 0; i < this->observed_values.size(); i++) {
471 distribution->observed_values[i] =
472 this->observed_values[i].initial_value_m;
473 }
474 // set relative info
475 distribution->expected_values.resize(this->expected_values.size());
476 for (size_t i = 0; i < this->expected_values.size(); i++) {
477 distribution->expected_values[i] =
478 this->expected_values[i].initial_value_m;
479 }
480 distribution->log_sd.resize(this->log_sd.size());
481 for (size_t i = 0; i < this->log_sd.size(); i++) {
482 distribution->log_sd[i] = this->log_sd[i].initial_value_m;
483 if (this->log_sd[i].estimation_type_m.get() == "fixed_effects") {
484 ss.str("");
485 ss << "dnorm." << this->id_m << ".log_sd." << this->log_sd[i].id_m;
486 info->RegisterParameterName(ss.str());
487 info->RegisterParameter(distribution->log_sd[i]);
488 }
489 if (this->log_sd[i].estimation_type_m.get() == "random_effects") {
490 FIMS_ERROR_LOG("standard deviations cannot be set to random effects");
491 }
492 }
493 info->variable_map[this->log_sd.id_m] = &(distribution)->log_sd;
494
495 distribution->use_mean = this->use_mean_m.get();
496 distribution->expected_mean.resize(this->expected_mean.size());
497 for (size_t i = 0; i < this->expected_mean.size(); i++) {
498 distribution->expected_mean[i] = this->expected_mean[i].initial_value_m;
499 if (this->expected_mean[i].estimation_type_m.get() == "fixed_effects") {
500 ss.str("");
501 ss << "dnorm." << this->id_m << ".expected_mean."
502 << this->expected_mean[i].id_m;
503 info->RegisterParameterName(ss.str());
504 info->RegisterParameter(distribution->expected_mean[i]);
505 }
506 if (this->expected_mean[i].estimation_type_m.get() == "random_effects") {
507 FIMS_ERROR_LOG("expected_mean cannot be set to random effects");
508 }
509 }
510 info->variable_map[this->expected_mean.id_m] =
512
513 info->density_components[distribution->id] = distribution;
514
515 return true;
516 }
517
522 virtual bool add_to_fims_tmb() {
525
526 return true;
527 }
528
529#endif
530};
531
537 public:
570
582
587
592 virtual uint32_t get_id() { return this->id_m; }
593
599 this->interface_observed_data_id_m.set(observed_data_id);
600 return true;
601 }
602
611 virtual bool set_distribution_links(std::string input_type,
612 Rcpp::IntegerVector ids) {
613 this->input_type_m.set(input_type);
614 this->key_m->resize(ids.size());
615 for (R_xlen_t i = 0; i < ids.size(); i++) {
616 this->key_m->at(i) = ids[i];
617 }
618 return true;
619 }
620
627 virtual double evaluate() {
629 dlnorm.observed_values.resize(this->observed_values.size());
630 dlnorm.expected_values.resize(this->expected_values.size());
631 dlnorm.log_sd.resize(this->log_sd.size());
632 // dlnorm.input_type = "prior";
633 for (size_t i = 0; i < this->observed_values.size(); i++) {
634 dlnorm.observed_values[i] = this->observed_values[i].initial_value_m;
635 }
636 for (size_t i = 0; i < this->expected_values.size(); i++) {
637 dlnorm.expected_values[i] = this->expected_values[i].initial_value_m;
638 }
639 for (size_t i = 0; i < this->log_sd.size(); i++) {
640 dlnorm.log_sd[i] = this->log_sd[i].initial_value_m;
641 }
642 return dlnorm.evaluate();
643 }
644
649 virtual void finalize() {
650 if (this->finalized) {
651 // log warning that finalize has been called more than once.
652 FIMS_WARNING_LOG("LogNormalLPDF " + fims::to_string(this->id_m) +
653 " has been finalized already.");
654 }
655
656 this->finalized = true; // indicate this has been called already
657
658 std::shared_ptr<fims_info::Information<double>> info =
660
662
663 // search for density component in Information
664 it = info->density_components.find(this->id_m);
665 // if not found, just return
666 if (it == info->density_components.end()) {
667 FIMS_WARNING_LOG("LogNormalLPDF " + fims::to_string(this->id_m) +
668 " not found in Information.");
669 return;
670 } else {
671 std::shared_ptr<fims_distributions::LogNormalLPDF<double>> dlnorm =
672 std::dynamic_pointer_cast<fims_distributions::LogNormalLPDF<double>>(
673 it->second);
674
675 this->lpdf_value = dlnorm->lpdf;
676
677 size_t n_x = dlnorm->get_n_x();
678
679 if (this->log_sd.size() != n_x) {
680 // If log_sd size == 1 (scalar), repeat the entry
681 if (this->log_sd.size() == 1) {
682 auto tmp = this->log_sd[0]; // copy the one log_sd param
683 this->log_sd.resize(n_x);
684 for (size_t i = 0; i < n_x; ++i) {
685 this->log_sd[i] = tmp; // copies all fields in Param
686 }
687 } else {
688 // Handle error
690 "log_sd size does not match number of observations and is not "
691 "scalar.");
692 }
693 }
694
695 for (size_t i = 0; i < n_x; i++) {
696 if (this->log_sd[i].estimation_type_m.get() == "constant") {
697 this->log_sd[i].final_value_m = this->log_sd[i].initial_value_m;
698 } else {
699 this->log_sd[i].final_value_m = dlnorm->log_sd.get_force_scalar(i);
700 }
701 }
702
703 this->lpdf_vec = RealVector(n_x);
704 if (this->expected_values.size() == 1) {
705 this->expected_values.resize(n_x);
706 }
707 if (this->observed_values.size() == 1) {
708 this->observed_values.resize(n_x);
709 }
710 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
711 this->lpdf_vec[i] = dlnorm->lpdf_vec[i];
712 this->expected_values[i].final_value_m = dlnorm->get_expected(i);
713 this->observed_values[i].final_value_m = dlnorm->get_observed(i);
714 }
715 }
716 }
717
725 virtual std::string to_json() {
726 std::stringstream ss;
727
728 ss << "{\n";
729 ss << " \"module_name\": \"density\",\n";
730 ss << " \"module_id\": " << this->id_m << ",\n";
731 ss << " \"module_type\": \"log_normal\",\n";
732 ss << " \"observed_data_id\" : " << this->interface_observed_data_id_m
733 << ",\n";
734 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
735 ss << " \"density_component\": {\n";
736 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
737 ss << " \"value\":[";
738 if (this->lpdf_vec.size() == 0) {
739 ss << "],\n";
740 } else {
741 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
742 ss << this->value_to_string(this->lpdf_vec[i]) << ", ";
743 }
744 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
745
746 ss << "],\n";
747 }
748 ss << " \"expected_values\":[";
749 if (this->expected_values.size() == 0) {
750 ss << "],\n";
751 } else {
752 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
753 ss << this->value_to_string(this->expected_values[i].final_value_m)
754 << ", ";
755 }
756 ss << this->value_to_string(
757 this->expected_values[this->expected_values.size() - 1]
758 .final_value_m);
759
760 ss << "],\n";
761 }
762 ss << " \"log_sd_values\":[";
763 if (this->log_sd.size() == 0) {
764 ss << "],\n";
765 } else {
766 for (R_xlen_t i = 0; i < static_cast<R_xlen_t>(this->log_sd.size()) - 1;
767 i++) {
768 ss << this->value_to_string(this->log_sd[i].final_value_m) << ", ";
769 }
770 ss << this->value_to_string(
771 this->log_sd[this->log_sd.size() - 1].final_value_m)
772 << "],\n";
773 }
774 ss << " \"observed_values\":[";
775 if (this->observed_values.size() == 0) {
776 ss << "]\n";
777 } else {
778 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
779 ss << this->observed_values[i].final_value_m << ", ";
780 }
781 ss << this->observed_values[this->observed_values.size() - 1]
782 .final_value_m
783 << "]\n";
784 }
785 ss << " }}\n";
786 return ss.str();
787 }
788
789#ifdef TMB_MODEL
790
791 template <typename Type>
793 std::shared_ptr<fims_info::Information<Type>> info =
795
796 std::shared_ptr<fims_distributions::LogNormalLPDF<Type>> distribution =
797 std::make_shared<fims_distributions::LogNormalLPDF<Type>>();
798
799 // set relative info
800 distribution->id = this->id_m;
801 std::stringstream ss;
802 distribution->observed_data_id_m = interface_observed_data_id_m;
803 distribution->input_type = this->input_type_m;
804 distribution->key.resize(this->key_m->size());
805 for (size_t i = 0; i < this->key_m->size(); i++) {
806 distribution->key[i] = this->key_m->at(i);
807 }
808 distribution->observed_values.resize(this->observed_values.size());
809 for (size_t i = 0; i < this->observed_values.size(); i++) {
810 distribution->observed_values[i] =
811 this->observed_values[i].initial_value_m;
812 }
813 // set relative info
814 distribution->expected_values.resize(this->expected_values.size());
815 for (size_t i = 0; i < this->expected_values.size(); i++) {
816 distribution->expected_values[i] =
817 this->expected_values[i].initial_value_m;
818 }
819 distribution->log_sd.resize(this->log_sd.size());
820 for (size_t i = 0; i < this->log_sd.size(); i++) {
821 distribution->log_sd[i] = this->log_sd[i].initial_value_m;
822 if (this->log_sd[i].estimation_type_m.get() == "fixed_effects") {
823 ss.str("");
824 ss << "dlnorm." << this->id_m << ".log_sd." << this->log_sd[i].id_m;
825 info->RegisterParameterName(ss.str());
826 info->RegisterParameter(distribution->log_sd[i]);
827 }
828 if (this->log_sd[i].estimation_type_m.get() == "random_effects") {
829 FIMS_ERROR_LOG("standard deviations cannot be set to random effects");
830 }
831 }
832 info->variable_map[this->log_sd.id_m] = &(distribution)->log_sd;
833
834 info->density_components[distribution->id] = distribution;
835
836 return true;
837 }
838
843 virtual bool add_to_fims_tmb() {
846
847 return true;
848 }
849
850#endif
851};
852
858 public:
885
895
908
917 virtual uint32_t get_id() { return this->id_m; }
918
924 this->interface_observed_data_id_m.set(observed_data_id);
925 return true;
926 }
927
936 virtual bool set_distribution_links(std::string input_type,
937 Rcpp::IntegerVector ids) {
938 this->input_type_m.set(input_type);
939 this->key_m->resize(ids.size());
940 for (R_xlen_t i = 0; i < ids.size(); i++) {
941 this->key_m->at(i) = ids[i];
942 }
943 return true;
944 }
945
951 void set_note(std::string note) { this->notes.set(note); }
952
958 virtual double evaluate() {
960 // Declare TMBVector in this scope
961 dmultinom.observed_values.resize(this->observed_values.size());
962 dmultinom.expected_values.resize(this->expected_values.size());
963 for (size_t i = 0; i < observed_values.size(); i++) {
964 dmultinom.observed_values[i] = this->observed_values[i].initial_value_m;
965 }
966 for (size_t i = 0; i < expected_values.size(); i++) {
967 dmultinom.expected_values[i] = this->expected_values[i].initial_value_m;
968 }
969 dmultinom.dims.resize(2);
970 dmultinom.dims[0] = this->dims[0];
971 dmultinom.dims[1] = this->dims[1];
972 return dmultinom.evaluate();
973 }
974
975 void finalize() {
976 if (this->finalized) {
977 // log warning that finalize has been called more than once.
978 FIMS_WARNING_LOG("DmultinomDistributions " +
979 fims::to_string(this->id_m) +
980 " has been finalized already.");
981 }
982
983 this->finalized = true; // indicate this has been called already
984
985 std::shared_ptr<fims_info::Information<double>> info =
987
989
990 // search for density component in Information
991 it = info->density_components.find(this->id_m);
992 // if not found, just return
993 if (it == info->density_components.end()) {
994 FIMS_WARNING_LOG("DmultinomDistributions " + fims::to_string(this->id_m) +
995 " not found in Information.");
996 return;
997 } else {
998 std::shared_ptr<fims_distributions::MultinomialLPMF<double>> dmultinom =
999 std::dynamic_pointer_cast<
1001
1002 this->lpdf_value = dmultinom->lpdf;
1003
1004 size_t n_x = dmultinom->lpdf_vec.size();
1005 this->lpdf_vec = Rcpp::NumericVector(n_x);
1006 if (this->expected_values.size() != n_x) {
1007 this->expected_values.resize(n_x);
1008 }
1009 if (this->observed_values.size() != n_x) {
1010 this->observed_values.resize(n_x);
1011 }
1012 for (size_t i = 0; i < this->lpdf_vec.size(); i++) {
1013 this->lpdf_vec[i] = dmultinom->lpdf_vec[i];
1014 this->expected_values[i].final_value_m = dmultinom->get_expected(i);
1015 if (dmultinom->input_type != "data") {
1016 this->observed_values[i].final_value_m = dmultinom->get_observed(i);
1017 }
1018 }
1019 if (dmultinom->input_type == "data") {
1020 dims.resize(2);
1021 dims[0] = dmultinom->dims[0];
1022 dims[1] = dmultinom->dims[1];
1023 for (size_t i = 0; i < dims[0]; i++) {
1024 for (size_t j = 0; j < dims[1]; j++) {
1025 size_t idx = (i * dims[1]) + j;
1026 this->observed_values[idx].final_value_m = dmultinom->get_observed(
1027 static_cast<size_t>(i), static_cast<size_t>(j));
1028 }
1029 }
1030 }
1031 }
1032 }
1033
1041 virtual std::string to_json() {
1042 std::stringstream ss;
1043
1044 ss << "{\n";
1045 ss << " \"module_name\": \"density\",\n";
1046 ss << " \"module_id\": " << this->id_m << ",\n";
1047 ss << " \"module_type\": \"multinomial\",\n";
1048 ss << "\"observed_data_id\" : " << this->interface_observed_data_id_m
1049 << ",\n";
1050 ss << " \"input_type\" : \"" << this->input_type_m << "\",\n";
1051 ss << " \"density_component\": {\n";
1052 ss << " \"lpdf_value\": " << this->lpdf_value << ",\n";
1053 ss << " \"value\":[";
1054 if (this->lpdf_vec.size() == 0) {
1055 ss << "],\n";
1056 } else {
1057 for (size_t i = 0; i < this->lpdf_vec.size() - 1; i++) {
1058 ss << this->value_to_string(this->lpdf_vec[i]);
1059 ss << ", ";
1060 }
1061 ss << this->value_to_string(this->lpdf_vec[this->lpdf_vec.size() - 1]);
1062
1063 ss << "],\n";
1064 }
1065 ss << " \"expected_values\":[";
1066 if (this->expected_values.size() == 0) {
1067 ss << "],\n";
1068 } else {
1069 for (size_t i = 0; i < this->expected_values.size() - 1; i++) {
1070 ss << this->value_to_string(this->expected_values[i].final_value_m)
1071 << ", ";
1072 }
1073 ss << this->value_to_string(
1074 this->expected_values[this->expected_values.size() - 1]
1075 .final_value_m);
1076
1077 ss << "],\n";
1078 }
1079 // no log_sd_values for multinomial
1080 ss << " \"observed_values\":[";
1081 if (this->observed_values.size() == 0) {
1082 ss << "]\n";
1083 } else {
1084 for (size_t i = 0; i < this->observed_values.size() - 1; i++) {
1085 ss << this->observed_values[i].final_value_m << ", ";
1086 }
1087 ss << this->observed_values[this->observed_values.size() - 1]
1088 .final_value_m
1089 << "]\n";
1090 }
1091 ss << " }}\n";
1092 return ss.str();
1093 }
1094
1095#ifdef TMB_MODEL
1096
1097 template <typename Type>
1099 std::shared_ptr<fims_info::Information<Type>> info =
1101
1102 std::shared_ptr<fims_distributions::MultinomialLPMF<Type>> distribution =
1103 std::make_shared<fims_distributions::MultinomialLPMF<Type>>();
1104
1105 distribution->id = this->id_m;
1106 distribution->observed_data_id_m = interface_observed_data_id_m;
1107 distribution->input_type = this->input_type_m;
1108 distribution->key.resize(this->key_m->size());
1109 for (size_t i = 0; i < this->key_m->size(); i++) {
1110 distribution->key[i] = this->key_m->at(i);
1111 }
1112 distribution->observed_values.resize(this->observed_values.size());
1113 for (size_t i = 0; i < this->observed_values.size(); i++) {
1114 distribution->observed_values[i] =
1115 this->observed_values[i].initial_value_m;
1116 }
1117 // set relative info
1118 distribution->expected_values.resize(this->expected_values.size());
1119 for (size_t i = 0; i < this->expected_values.size(); i++) {
1120 distribution->expected_values[i] =
1121 this->expected_values[i].initial_value_m;
1122 }
1123
1124 info->density_components[distribution->id] = distribution;
1125 return true;
1126 }
1127
1128 virtual bool add_to_fims_tmb() {
1131
1132 return true;
1133 }
1134
1135#endif
1136};
1137
1138#endif
Rcpp interface that serves as the parent class for Rcpp distribution interfaces. This type should be ...
Definition rcpp_distribution.hpp:19
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:154
double lpdf_value
The log probability density function value.
Definition rcpp_distribution.hpp:76
virtual ~DistributionsInterfaceBase()
The destructor.
Definition rcpp_distribution.hpp:103
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:32
DistributionsInterfaceBase(const DistributionsInterfaceBase &other)
Construct a new Distributions Interface Base object.
Definition rcpp_distribution.hpp:93
uint32_t id_m
The local ID of the DistributionsInterfaceBase object.
Definition rcpp_distribution.hpp:28
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:60
SharedInt interface_observed_data_id_m
The ID of the observed data object, which is set to -999.
Definition rcpp_distribution.hpp:71
virtual bool set_distribution_mean(double input_value)
Set the expected mean value for the distribution.
Definition rcpp_distribution.hpp:146
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:118
static uint32_t id_g
The static ID of the DistributionsInterfaceBase object.
Definition rcpp_distribution.hpp:24
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:67
SharedString input_type_m
The type of density input. The options are prior, re, or data.
Definition rcpp_distribution.hpp:36
DistributionsInterfaceBase()
The constructor.
Definition rcpp_distribution.hpp:80
The Rcpp interface for Dlnorm to instantiate from R: dlnorm_ <- methods::new(DlnormDistribution).
Definition rcpp_distribution.hpp:536
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:611
DlnormDistributionsInterface(const DlnormDistributionsInterface &other)
Construct a new Dlnorm Distributions Interface object.
Definition rcpp_distribution.hpp:576
DlnormDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:564
ParameterVector expected_values
The expected values, which would be the mean of log(x) for this distribution.
Definition rcpp_distribution.hpp:546
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:725
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:592
virtual double evaluate()
Evaluate lognormal probability density function (pdf). The natural log of the pdf is returned.
Definition rcpp_distribution.hpp:627
ParameterVector observed_values
Observed data.
Definition rcpp_distribution.hpp:541
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:554
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:598
virtual void finalize()
Extracts the derived quantities from Information to the Rcpp object.
Definition rcpp_distribution.hpp:649
virtual ~DlnormDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:586
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:559
The Rcpp interface for Dmultinom to instantiate from R: dmultinom_ <- methods::new(DmultinomDistribut...
Definition rcpp_distribution.hpp:857
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:936
SharedString notes
TODO: document this.
Definition rcpp_distribution.hpp:884
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:917
ParameterVector expected_values
The expected values, which should be a vector of length K where each value specifies the probability ...
Definition rcpp_distribution.hpp:868
void finalize()
Extracts derived quantities back to the Rcpp interface object from the Information object.
Definition rcpp_distribution.hpp:975
virtual double evaluate()
Definition rcpp_distribution.hpp:958
ParameterVector observed_values
Observed data, which should be a vector of length K of integers.
Definition rcpp_distribution.hpp:862
void set_note(std::string note)
Set the note object.
Definition rcpp_distribution.hpp:951
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:923
RealVector dims
The dimensions of the number of rows and columns of the multivariate dataset.
Definition rcpp_distribution.hpp:873
DmultinomDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:889
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:1041
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:878
DmultinomDistributionsInterface(const DmultinomDistributionsInterface &other)
Construct a new Dmultinom Distributions Interface object.
Definition rcpp_distribution.hpp:901
virtual ~DmultinomDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:912
The Rcpp interface for Dnorm to instantiate from R: dnorm_ <- methods::new(DnormDistribution).
Definition rcpp_distribution.hpp:173
DnormDistributionsInterface(const DnormDistributionsInterface &other)
Construct a new Dnorm Distributions Interface object.
Definition rcpp_distribution.hpp:215
ParameterVector expected_mean
The expected mean, which would be the mean of x for this distribution.
Definition rcpp_distribution.hpp:188
DnormDistributionsInterface()
The constructor.
Definition rcpp_distribution.hpp:203
virtual std::string to_json()
Converts the data to json representation for the output.
Definition rcpp_distribution.hpp:385
virtual ~DnormDistributionsInterface()
The destructor.
Definition rcpp_distribution.hpp:226
RealVector lpdf_vec
Vector that records the individual log probability function for each observation.
Definition rcpp_distribution.hpp:198
virtual double evaluate()
Evaluate normal probability density function (pdf). The natural log of the pdf is returned.
Definition rcpp_distribution.hpp:272
ParameterVector expected_values
The expected values, which would be the mean of x for this distribution.
Definition rcpp_distribution.hpp:183
virtual void finalize()
Extracts the derived quantities from Information to the Rcpp object.
Definition rcpp_distribution.hpp:298
virtual bool set_observed_data(int observed_data_id)
Set the unique ID for the observed data object.
Definition rcpp_distribution.hpp:238
ParameterVector observed_values
Observed data.
Definition rcpp_distribution.hpp:178
virtual bool set_distribution_mean(double input_value)
Set the expected mean value for the distribution.
Definition rcpp_distribution.hpp:246
ParameterVector log_sd
The uncertainty, which would be the standard deviation of x for the normal distribution.
Definition rcpp_distribution.hpp:193
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:256
virtual uint32_t get_id()
Gets the ID of the interface base object.
Definition rcpp_distribution.hpp:232
Base class for all interface objects.
Definition rcpp_interface_base.hpp:575
bool finalized
Is the object already finalized? The default is false.
Definition rcpp_interface_base.hpp:580
std::string value_to_string(double value)
Report the parameter value as a string.
Definition rcpp_interface_base.hpp:616
static std::vector< std::shared_ptr< FIMSRcppInterfaceBase > > fims_interface_objects
FIMS interface object vectors.
Definition rcpp_interface_base.hpp:585
virtual bool add_to_fims_tmb()
A virtual method to inherit to add objects to the TMB model.
Definition rcpp_interface_base.hpp:590
An Rcpp interface class that defines the ParameterVector class.
Definition rcpp_interface_base.hpp:143
void resize(size_t size)
Resizes a ParameterVector to the desired length.
Definition rcpp_interface_base.hpp:288
size_t size()
Returns the size of a ParameterVector.
Definition rcpp_interface_base.hpp:281
uint32_t id_m
The local ID of the Parameter object.
Definition rcpp_interface_base.hpp:156
Parameter & get(size_t pos)
An internal accessor for calling a position of a ParameterVector from R.
Definition rcpp_interface_base.hpp:260
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:276
An Rcpp interface class that defines the RealVector class.
Definition rcpp_interface_base.hpp:375
size_t size()
Returns the size of a RealVector.
Definition rcpp_interface_base.hpp:536
void resize(size_t size)
Resizes a RealVector to the desired length.
Definition rcpp_interface_base.hpp:543
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:639
#define FIMS_ERROR_LOG(MESSAGE)
Definition def.hpp:655
void clear_internal()
Clears the internal objects.
Definition rcpp_interface.hpp:235
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:32
virtual const Type evaluate()
Evaluates the normal log probability density function.
Definition normal_lpdf.hpp:59