deal.II version GIT relicensing-2173-gae8fc9d14b 2024-11-24 06:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
MCMC-Laplace.h
Go to the documentation of this file.
1 5,
486 /* fe_degree = */ 1,
487 "");
488 LogLikelihood::Gaussian log_likelihood(exact_solution, 0.05);
489 LogPrior::LogGaussian log_prior(0, 2);
490
491 const unsigned int n_theta = 64;
492 for (unsigned int test=0; test<10; ++test)
493 {
494 std::cout << "Generating output for test " << test << std::endl;
495
496 // For each test, read the input file...
497 std::ifstream test_input ("../testing/input." + std::to_string(test) + ".txt");
498 Assert (test_input, ExcIO());
499
500 Vector<double> theta(n_theta);
501 for (unsigned int i=0; i<n_theta; ++i)
502 test_input >> theta[i];
503
504 // ...run it through the forward simulator to get the
505 // z vector (which we then output to a file)...
506 const Vector<double> z = laplace_problem.evaluate(theta);
507 std::ofstream test_output_z ("output." + std::to_string(test) + ".z.txt");
508 z.print(test_output_z, 16);
509
510 // ...and then also evaluate prior and likelihood for these
511 // theta vectors:
512 std::ofstream test_output_likelihood ("output." + std::to_string(test) + ".loglikelihood.txt");
513 test_output_likelihood.precision(12);
514 test_output_likelihood << log_likelihood.log_likelihood(z) << std::endl;
515
516 std::ofstream test_output_prior ("output." + std::to_string(test) + ".logprior.txt");
517 test_output_prior.precision(12);
518 test_output_prior << log_prior.log_prior(theta) << std::endl;
519 }
520}
521@endcode
522This code reads in each of the input files (assuming that the executable is located in a
523build directory parallel to the `testing/` directory) and outputs its results into the
524current directory. The inputs you get from a modified program should then be compared
525against the ones stored in the `testing/` directory. They should match to several digits.
526
527
528
529An alternative implementation in Matlab
530---------------------------------------
531
532To facilitate experiments, this directory also contains alternative
533implementations of the benchmark. The first one was written by David
534Aristoff in Matlab and can be found in the `Matlab/` directory. As is
535common in Matlab programs, each function is provided in its own
536file. We have verified that the program generates the same results (to
53712 or more digits) as the C++ program, using the tests mentioned in
538the previous section.
539
540
541
542An alternative implementation in Python
543---------------------------------------
544
545Another alternative, written by Wolfgang Bangerth, can be found in the
546`Python/` directory and uses Python3. As is common for Python
547programs, the whole program is provided in one file. As for the Matlab
548version, we have verified that the program generates the same results
549(to 12 or more digits) as the C++ program, using the tests mentioned
550in the previous section. In fact, if you just execute the program
551as-is, it runs diagnostics that output the errors.
552
553This Python version is essentially a literal translation of the Matlab
554code. It is not nearly as efficient (taking around 5 times as much
555time for each function evaluation) and could probably be optimized
556substantially if desired. A good starting point is the insertion of
557the local elements into the global matrix in the line
558@code{.sh}
559 A[np.ix_(dof,dof)] += theta_loc * A_loc
560@endcode
561that takes up a substantial fraction of the overall run time.
562
563
564<a name="ann-Matlab/exact_values.m"></a>
565<h1>Annotated version of Matlab/exact_values.m</h1>
566@code{.m}
567%% -----------------------------------------------------------------------------
568%%
569%% SPDX-License-Identifier: LGPL-2.1-or-later
570%% Copyright (C) 2022 by Wolfgang Bangerth
571%%
572%% This file is part of the deal.II code gallery.
573%%
574%% -----------------------------------------------------------------------------
575
576%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
577%%%%%%%%%%%% list of "exact" measurement values, z_hat %%%%%%%%%%%%%%%%%%%%
578%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
579
580z_hat = [0.06076511762259369;
581 0.09601910120848481;
582 0.1238852517838584;
583 0.1495184117375201;
584 0.1841596127549784;
585 0.2174525028261122;
586 0.2250996160898698;
587 0.2197954769002993;
588 0.2074695698370926;
589 0.1889996477663016;
590 0.1632722532153726;
591 0.1276782480038186;
592 0.07711845915789312;
593 0.09601910120848552;
594 0.2000589533367983;
595 0.3385592591951766;
596 0.3934300024647806;
597 0.4040223892461541;
598 0.4122329537843092;
599 0.4100480091545554;
600 0.3949151637189968;
601 0.3697873264791232;
602 0.33401826235924;
603 0.2850397806663382;
604 0.2184260032478671;
605 0.1271121156350957;
606 0.1238852517838611;
607 0.3385592591951819;
608 0.7119285162766475;
609 0.8175712861756428;
610 0.6836254116578105;
611 0.5779452419831157;
612 0.5555615956136897;
613 0.5285181561736719;
614 0.491439702849224;
615 0.4409367494853282;
616 0.3730060082060772;
617 0.2821694983395214;
618 0.1610176733857739;
619 0.1495184117375257;
620 0.3934300024647929;
621 0.8175712861756562;
622 0.9439154625527653;
623 0.8015904115095128;
624 0.6859683749254024;
625 0.6561235366960599;
626 0.6213197201867315;
627 0.5753611315000049;
628 0.5140091754526823;
629 0.4325325506354165;
630 0.3248315148915482;
631 0.1834600412730086;
632 0.1841596127549917;
633 0.4040223892461832;
634 0.6836254116578439;
635 0.8015904115095396;
636 0.7870119561144977;
637 0.7373108331395808;
638 0.7116558878070463;
639 0.6745179049094283;
640 0.6235300574156917;
641 0.5559332704045935;
642 0.4670304994474178;
643 0.3499809143811;
644 0.19688263746294;
645 0.2174525028261253;
646 0.4122329537843404;
647 0.5779452419831566;
648 0.6859683749254372;
649 0.7373108331396063;
650 0.7458811983178246;
651 0.7278968022406559;
652 0.6904793535357751;
653 0.6369176452710288;
654 0.5677443693743215;
655 0.4784738764865867;
656 0.3602190632823262;
657 0.2031792054737325;
658 0.2250996160898818;
659 0.4100480091545787;
660 0.5555615956137137;
661 0.6561235366960938;
662 0.7116558878070715;
663 0.727896802240657;
664 0.7121928678670187;
665 0.6712187391428729;
666 0.6139157775591492;
667 0.5478251665295381;
668 0.4677122687599031;
669 0.3587654911000848;
670 0.2050734291675918;
671 0.2197954769003094;
672 0.3949151637190157;
673 0.5285181561736911;
674 0.6213197201867471;
675 0.6745179049094407;
676 0.690479353535786;
677 0.6712187391428787;
678 0.6178408289359514;
679 0.5453605027237883;
680 0.489575966490909;
681 0.4341716881061278;
682 0.3534389974779456;
683 0.2083227496961347;
684 0.207469569837099;
685 0.3697873264791366;
686 0.4914397028492412;
687 0.5753611315000203;
688 0.6235300574157017;
689 0.6369176452710497;
690 0.6139157775591579;
691 0.5453605027237935;
692 0.4336604929612851;
693 0.4109641743019312;
694 0.3881864790111245;
695 0.3642640090182592;
696 0.2179599909280145;
697 0.1889996477663011;
698 0.3340182623592461;
699 0.4409367494853381;
700 0.5140091754526943;
701 0.5559332704045969;
702 0.5677443693743304;
703 0.5478251665295453;
704 0.4895759664908982;
705 0.4109641743019171;
706 0.395727260284338;
707 0.3778949322004734;
708 0.3596268271857124;
709 0.2191250268948948;
710 0.1632722532153683;
711 0.2850397806663325;
712 0.373006008206081;
713 0.4325325506354207;
714 0.4670304994474315;
715 0.4784738764866023;
716 0.4677122687599041;
717 0.4341716881061055;
718 0.388186479011099;
719 0.3778949322004602;
720 0.3633362567187364;
721 0.3464457261905399;
722 0.2096362321365655;
723 0.1276782480038148;
724 0.2184260032478634;
725 0.2821694983395252;
726 0.3248315148915535;
727 0.3499809143811097;
728 0.3602190632823333;
729 0.3587654911000799;
730 0.3534389974779268;
731 0.3642640090182283;
732 0.35962682718569;
733 0.3464457261905295;
734 0.3260728953424643;
735 0.180670595355394;
736 0.07711845915789244;
737 0.1271121156350963;
738 0.1610176733857757;
739 0.1834600412730144;
740 0.1968826374629443;
741 0.2031792054737354;
742 0.2050734291675885;
743 0.2083227496961245;
744 0.2179599909279998;
745 0.2191250268948822;
746 0.2096362321365551;
747 0.1806705953553887;
748 0.1067965550010013];
749@endcode
750
751
752<a name="ann-Matlab/forward_solver.m"></a>
753<h1>Annotated version of Matlab/forward_solver.m</h1>
754@code{.m}
755%% -----------------------------------------------------------------------------
756%%
757%% SPDX-License-Identifier: LGPL-2.1-or-later
758%% Copyright (C) 2022 by Wolfgang Bangerth
759%%
760%% This file is part of the deal.II code gallery.
761%%
762%% -----------------------------------------------------------------------------
763
764%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
765%%%%%%%%%%%%%%%%%%%%%% forward solver function %%%%%%%%%%%%%%%%%%%%%%%%%%%%
766%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
767%INPUTS:
768%theta = current 8x8 parameter matrix
769%lbl = cell labeling function
770%A_loc = matrix of local contributions to A
771%Id = Identity matrix of size 128x128
772%boundaries = labels of boundary cells
773%b = right hand side of linear system (AU = b)
774%M = measurement matrix
775%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
776%OUTPUTS:
777%z = vector of measurements
778%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
779
780function z = forward_solver(theta,lbl,A_loc,Id,boundaries,b,M)
781
782%initialize matrix A for FEM linear solve, AU = b
783A = zeros(33^2,33^2);
784
785%loop over cells to build A
786for i=0:31
787 for j=0:31 %build A by summing over contribution from each cell
788
789 %find local coefficient in 8x8 grid
790 theta_loc = theta(floor(i/4)+1,floor(j/4)+1);
791
792 %update A by including contribution from cell (i,j)
793 dof = [lbl(i,j),lbl(i,j+1),lbl(i+1,j+1),lbl(i+1,j)];
794 A(dof,dof) = A(dof,dof) + theta_loc*A_loc;
795 end
796end
797
798%enforce boundary condition
799A(boundaries,:) = Id(boundaries,:);
800A(:,boundaries) = Id(:,boundaries);
801
802%sparsify A
803A = sparse(A);
804
805%solve linear equation for coefficients, U
806U = A\b;
807
808%get new z values
809z = M*U;
810@endcode
811
812
813<a name="ann-Matlab/get_statistics.m"></a>
814<h1>Annotated version of Matlab/get_statistics.m</h1>
815@code{.m}
816%% -----------------------------------------------------------------------------
817%%
818%% SPDX-License-Identifier: LGPL-2.1-or-later
819%% Copyright (C) 2022 by Wolfgang Bangerth
820%%
821%% This file is part of the deal.II code gallery.
822%%
823%% -----------------------------------------------------------------------------
824
825%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
826%%%%%%%%%%%%%%%%% compute statistics on data set %%%%%%%%%%%%%%%%%%%%%%%%%%
827%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
828%INPUTS:
829%data = tensor of theta samples from each lag time and chain
830%theta_means = means of theta over each independent chain
831%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
832%OUTPUTS:
833%theta_mean = overall mean of chains
834%covars = covariance matrices of each independent chain
835%autocovar = mean of autocovariance matrix over all the chains
836%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
837
838function [theta_mean,covars,autocovar] = get_statistics(data,theta_means);
839
840%compute overall mean of data, and get size of data matrix
841theta_mean = mean(theta_means,3);
842[~,~,L,N] = size(data);
843
844%initialize covariance matrices and mean autocovariance matrix
845covars = zeros(64,64,N);
846autocovar = zeros(64,64,2*L-1);
847
848%compute covariance matrices and mean autocovariance matrix
849for n=1:N %loop over independent Markov chains
850
851 %get data from chain n
852 data_ = reshape(permute(data(:,:,:,n),[3 2 1]),[L 64]);
853
854 %compute autocovariance matrix of chain n
855 mat = xcov(data_,'unbiased');
856
857 %store covariance matrix of chain n
858 covars(:,:,n) = reshape(mat(L,:),64,64);
859
860 %update mean autocovariance matrix
861 autocovar = autocovar + reshape(mat',[64 64 2*L-1]);
862end
863
864%compute mean of autocovariance matrix
865autocovar = autocovar(1:64,1:64,L:2*L-1)/N;
866@endcode
867
868
869<a name="ann-Matlab/log_probability.m"></a>
870<h1>Annotated version of Matlab/log_probability.m</h1>
871@code{.m}
872%% -----------------------------------------------------------------------------
873%%
874%% SPDX-License-Identifier: LGPL-2.1-or-later
875%% Copyright (C) 2022 by Wolfgang Bangerth
876%%
877%% This file is part of the deal.II code gallery.
878%%
879%% -----------------------------------------------------------------------------
880
881%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
882%%%%%%%%%%%%%%%%% compute log probability, log pi %%%%%%%%%%%%%%%%%%%%%%%%%
883%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
884%INPUTS:
885%theta = current 8x8 parameter matrix
886%z = current vector of measurements
887%z_hat = vector of "exact" measurements
888%sig = standard deviation parameter in likelihood
889%sig_pr = standard deviation parameter in prior
890%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
891%OUTPUTS:
892%log_pi = logarithm of posterior probability
893%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
894
895function log_pi = log_probability(theta,z,z_hat,sig,sig_pr)
896
897%compute log likelihood
898log_L = -sum((z-z_hat).^2)/(2*sig^2);
899
900%compute log prior
901log_pi_pr = -sum(log(theta).^2,'all')/(2*sig_pr^2);
902
903%compute log posterior
904log_pi = log_L+log_pi_pr;
905@endcode
906
907
908<a name="ann-Matlab/main.m"></a>
909<h1>Annotated version of Matlab/main.m</h1>
910@code{.m}
911%% -----------------------------------------------------------------------------
912%%
913%% SPDX-License-Identifier: LGPL-2.1-or-later
914%% Copyright (C) 2022 by Wolfgang Bangerth
915%%
916%% This file is part of the deal.II code gallery.
917%%
918%% -----------------------------------------------------------------------------
919
920%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
921%%%%%%%%% Run MCMC sampler to estimate posterior distribution %%%%%%%%%%%%%
922%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
923
924%define number of chains, chain length, and lag time
925N = input('number of independent Markov chains: ');
926N_L = input('length of each Markov chain: ');
927lag = input('lag time for measurements: ');
928workers = input('number of parallel workers: ');
929L = N_L/lag;
930
931%open Matlab parallel pool
932parpool(workers)
933
934%load precomputations
935load precomputations.mat
936
937%define lag time and data matrix
938data = zeros(8,8,L,N); %data matrix of samples at lag times
939theta_means = zeros(8,8,N); %overall mean of theta
940
941tic
942
943parfor n=1:N
944
945 %set initial theta, theta mean, and z values of chain
946 theta = theta0;
947 theta_mean = 0;
948 z = z0;
949
950 for m=1:L
951
952 for l=1:lag
953
954 %define proposal, theta_tilde
955 xi = normrnd(0,sig_prop,[8 8]);
956 theta_tilde = theta.*exp(xi);
957
958 %compute new z values
959 z_tilde = forward_solver_(theta_tilde);
960
961 %compute posterior log probability of theta_tilde
962 log_pi_tilde = log_probability_(theta_tilde,z_tilde);
963 log_pi = log_probability_(theta,z);
964
965 %compute acceptance probability; accept proposal appropriately
966 accept = exp(log_pi_tilde-log_pi) ...
967 *prod(theta_tilde./theta,'all');
968 if rand<accept
969 theta = theta_tilde; %accept new theta values
970 z = z_tilde; %record associated measurements
971 end
972
973 %update mean of theta
974 theta_mean = theta_mean + theta;
975
976 end
977
978 %update data matrix
979 data(:,:,m,n) = theta;
980
981 end
982
983 %update theta means
984 theta_means(:,:,n) = theta_mean/N_L;
985
986end
987
988toc
989
990%compute statistics on data set
991[theta_mean,covars,autocovar] = get_statistics(data,theta_means);
992
993%save data to Matlab workspace, labeled by N and N_L
994save (['data_N_' num2str(N) '_N_L_ ' num2str(N_L) '.mat'])
995@endcode
996
997
998<a name="ann-Matlab/plot_solution.m"></a>
999<h1>Annotated version of Matlab/plot_solution.m</h1>
1000@code{.m}
1001%% -----------------------------------------------------------------------------
1002%%
1003%% SPDX-License-Identifier: LGPL-2.1-or-later
1004%% Copyright (C) 2022 by Wolfgang Bangerth
1005%%
1006%% This file is part of the deal.II code gallery.
1007%%
1008%% -----------------------------------------------------------------------------
1009
1010%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1011%%%%%%%%%%%%%%% plots solution, u, to Poisson equation %%%%%%%%%%%%%%%%%%%%
1012%%%%%%%%%%%%%%%% associated to theta and a 32x32 mesh %%%%%%%%%%%%%%%%%%%%%
1013%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1014
1015%input matrix theta for plotting (e.g., theta = theta_hat)
1016theta = input('choose theta for plotting u: theta = ');
1017
1018%construct mass matrix, M_plot, for plotting
1019xsp = 0:0.02:1;
1020n = length(xsp);
1021Mp = zeros(n,n,33^2);
1022for k=1:33^2
1023 c = inv_lbl(k);
1024 for i=1:n
1025 for j=1:n
1026 Mp(i,j,k) = phi(xsp(i)-h*c(1),xsp(j)-h*c(2));
1027 end
1028 end
1029end
1030Mp = reshape(Mp,[n^2 33^2]);
1031
1032%run forward solver on mean of theta
1033A = zeros(33^2,33^2);
1034for i=0:31
1035 for j=0:31 %build A by summing over contribution from each cell
1036
1037 %find local coefficient in 8x8 grid
1038 theta_loc = theta(floor(i/4)+1,floor(j/4)+1);
1039
1040 %update A by including contribution from cell (i,j)
1041 dof = [lbl(i,j),lbl(i,j+1),lbl(i+1,j+1),lbl(i+1,j)];
1042 A(dof,dof) = A(dof,dof) + theta_loc*A_loc;
1043 end
1044end
1045
1046%enforce boundary condition
1047A(boundaries,:) = Id(boundaries,:);
1048
1049%sparsify A
1050A = sparse(A);
1051
1052%solve linear equation for coefficients, U
1053U = A\b;
1054
1055%close all current plots
1056close all
1057
1058%plot solution
1059figure
1060zs = reshape(Mp*U,[n n]);
1061surf(xsp,xsp,zs)
1062@endcode
1063
1064
1065<a name="ann-Matlab/precomputations.m"></a>
1066<h1>Annotated version of Matlab/precomputations.m</h1>
1067@code{.m}
1068%% -----------------------------------------------------------------------------
1069%%
1070%% SPDX-License-Identifier: LGPL-2.1-or-later
1071%% Copyright (C) 2022 by Wolfgang Bangerth
1072%%
1073%% This file is part of the deal.II code gallery.
1074%%
1075%% -----------------------------------------------------------------------------
1076
1077%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1078%%%%%%% do all precomputations necessary for MCMC simulations %%%%%%%%%%%%%
1079%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1080
1081%define mesh width
1082h = 1/32;
1083
1084%define characteristic function of unit square
1085S = @(x,y) heaviside(x).*heaviside(y) ...
1086 .*(1-heaviside(x-h)).*(1-heaviside(y-h));
1087
1088%define tent function on the domain [-h,h]x[-h,h]
1089phi = @(x,y) ((x+h).*(y+h).*S(x+h,y+h) + (h-x).*(h-y).*S(x,y) ...
1090 + (x+h).*(h-y).*S(x+h,y) + (h-x).*(y+h).*S(x,y+h))/h^2;
1091
1092%define function that converts from (i,j) to dof, and its inverse
1093lbl = @(i,j) 33*j+i+1;
1094inv_lbl = @(k) [k-1-33*floor((k-1)/33),floor((k-1)/33)];
1095
1096%construct measurement matrix, M
1097xs = 1/14:1/14:13/14; %measurement points
1098M = zeros(13,13,33^2);
1099for k=1:33^2
1100 c = inv_lbl(k);
1101 for i=1:13
1102 for j=1:13
1103 M(i,j,k) = phi(xs(i)-h*c(1),xs(j)-h*c(2));
1104 end
1105 end
1106end
1107M = reshape(M,[13^2 33^2]);
1108
1109%construct exact coefficient matrix, theta_hat
1110theta_hat = ones(8,8);
1111theta_hat(2:3,2:3) = 0.1;
1112theta_hat(6:7,6:7) = 10;
1113
1114%construct local overlap matrix, A_loc, and identity matrix Id
1115A_loc = [2/3 -1/6 -1/3 -1/6;
1116 -1/6 2/3 -1/6 -1/3;
1117 -1/3 -1/6 2/3 -1/6;
1118 -1/6 -1/3 -1/6 2/3];
1119Id = eye(33^2);
1120
1121%locate boundary labels
1122boundaries = [lbl(0:1:32,0),lbl(0:1:32,32),lbl(0,1:1:31),lbl(32,1:1:31)];
1123
1124%define RHS of FEM linear system, AU = b
1125b = ones(33^2,1)*10*h^2;
1126b(boundaries) = zeros(128,1); %enforce boundary conditions on b
1127
1128%load exact z_hat values
1129exact_values
1130
1131%set global parameters and functions for simulation
1132sig = 0.05; %likelihood standard deviation
1133sig_pr = 2; %prior (log) standard deviation
1134sig_prop = 0.0725; %proposal (log) standard deviation
1135theta0 = ones(8,8); %initial theta values
1136forward_solver_ = @(theta) ...
1137 forward_solver(theta,lbl,A_loc,Id,boundaries,b,M);
1138log_probability_ = @(theta,z) log_probability(theta,z,z_hat,sig,sig_pr);
1139
1140%find initial z values
1141z0 = forward_solver_(theta0);
1142
1143save precomputations.mat
1144@endcode
1145
1146
1147<a name="ann-Python/forward_solver.py"></a>
1148<h1>Annotated version of Python/forward_solver.py</h1>
1149@code{.py}
1150## -----------------------------------------------------------------------------
1151##
1152## SPDX-License-Identifier: LGPL-2.1-or-later
1153## Copyright (C) 2022 by Wolfgang Bangerth
1154##
1155## This file is part of the deal.II code gallery.
1156##
1157## -----------------------------------------------------------------------------
1158
1159import numpy as np
1160import scipy.sparse
1161from scipy.sparse.linalg import spsolve
1162import time
1163
1164###########################################################################
1165############ list of "exact" measurement values, z_hat ####################
1166###########################################################################
1167
1168z_hat = np.array(
1169 [0.06076511762259369, 0.09601910120848481,
1170 0.1238852517838584, 0.1495184117375201,
1171 0.1841596127549784, 0.2174525028261122,
1172 0.2250996160898698, 0.2197954769002993,
1173 0.2074695698370926, 0.1889996477663016,
1174 0.1632722532153726, 0.1276782480038186,
1175 0.07711845915789312, 0.09601910120848552,
1176 0.2000589533367983, 0.3385592591951766,
1177 0.3934300024647806, 0.4040223892461541,
1178 0.4122329537843092, 0.4100480091545554,
1179 0.3949151637189968, 0.3697873264791232,
1180 0.33401826235924, 0.2850397806663382,
1181 0.2184260032478671, 0.1271121156350957,
1182 0.1238852517838611, 0.3385592591951819,
1183 0.7119285162766475, 0.8175712861756428,
1184 0.6836254116578105, 0.5779452419831157,
1185 0.5555615956136897, 0.5285181561736719,
1186 0.491439702849224, 0.4409367494853282,
1187 0.3730060082060772, 0.2821694983395214,
1188 0.1610176733857739, 0.1495184117375257,
1189 0.3934300024647929, 0.8175712861756562,
1190 0.9439154625527653, 0.8015904115095128,
1191 0.6859683749254024, 0.6561235366960599,
1192 0.6213197201867315, 0.5753611315000049,
1193 0.5140091754526823, 0.4325325506354165,
1194 0.3248315148915482, 0.1834600412730086,
1195 0.1841596127549917, 0.4040223892461832,
1196 0.6836254116578439, 0.8015904115095396,
1197 0.7870119561144977, 0.7373108331395808,
1198 0.7116558878070463, 0.6745179049094283,
1199 0.6235300574156917, 0.5559332704045935,
1200 0.4670304994474178, 0.3499809143811,
1201 0.19688263746294, 0.2174525028261253,
1202 0.4122329537843404, 0.5779452419831566,
1203 0.6859683749254372, 0.7373108331396063,
1204 0.7458811983178246, 0.7278968022406559,
1205 0.6904793535357751, 0.6369176452710288,
1206 0.5677443693743215, 0.4784738764865867,
1207 0.3602190632823262, 0.2031792054737325,
1208 0.2250996160898818, 0.4100480091545787,
1209 0.5555615956137137, 0.6561235366960938,
1210 0.7116558878070715, 0.727896802240657,
1211 0.7121928678670187, 0.6712187391428729,
1212 0.6139157775591492, 0.5478251665295381,
1213 0.4677122687599031, 0.3587654911000848,
1214 0.2050734291675918, 0.2197954769003094,
1215 0.3949151637190157, 0.5285181561736911,
1216 0.6213197201867471, 0.6745179049094407,
1217 0.690479353535786, 0.6712187391428787,
1218 0.6178408289359514, 0.5453605027237883,
1219 0.489575966490909, 0.4341716881061278,
1220 0.3534389974779456, 0.2083227496961347,
1221 0.207469569837099, 0.3697873264791366,
1222 0.4914397028492412, 0.5753611315000203,
1223 0.6235300574157017, 0.6369176452710497,
1224 0.6139157775591579, 0.5453605027237935,
1225 0.4336604929612851, 0.4109641743019312,
1226 0.3881864790111245, 0.3642640090182592,
1227 0.2179599909280145, 0.1889996477663011,
1228 0.3340182623592461, 0.4409367494853381,
1229 0.5140091754526943, 0.5559332704045969,
1230 0.5677443693743304, 0.5478251665295453,
1231 0.4895759664908982, 0.4109641743019171,
1232 0.395727260284338, 0.3778949322004734,
1233 0.3596268271857124, 0.2191250268948948,
1234 0.1632722532153683, 0.2850397806663325,
1235 0.373006008206081, 0.4325325506354207,
1236 0.4670304994474315, 0.4784738764866023,
1237 0.4677122687599041, 0.4341716881061055,
1238 0.388186479011099, 0.3778949322004602,
1239 0.3633362567187364, 0.3464457261905399,
1240 0.2096362321365655, 0.1276782480038148,
1241 0.2184260032478634, 0.2821694983395252,
1242 0.3248315148915535, 0.3499809143811097,
1243 0.3602190632823333, 0.3587654911000799,
1244 0.3534389974779268, 0.3642640090182283,
1245 0.35962682718569, 0.3464457261905295,
1246 0.3260728953424643, 0.180670595355394,
1247 0.07711845915789244, 0.1271121156350963,
1248 0.1610176733857757, 0.1834600412730144,
1249 0.1968826374629443, 0.2031792054737354,
1250 0.2050734291675885, 0.2083227496961245,
1251 0.2179599909279998, 0.2191250268948822,
1252 0.2096362321365551, 0.1806705953553887,
1253 0.1067965550010013])
1254
1255
1256###########################################################################
1257####### do all precomputations necessary for MCMC simulations #############
1258###########################################################################
1259
1260# Define the mesh width
1261h = 1/32
1262
1263# Define characteristic function of unit square
1264def heaviside(x) :
1265 if x<0 :
1266 return 0
1267 else :
1268 return 1
1269
1270def S(x,y) :
1271 return heaviside(x)*heaviside(y) * (1-heaviside(x-h))*(1-heaviside(y-h));
1272
1273# Define tent function on the domain [0,2h]x[0,2h]
1274def phi(x,y) :
1275 return ((x+h)*(y+h)*S(x+h,y+h) + (h-x)*(h-y)*S(x,y)
1276 + (x+h)*(h-y)*S(x+h,y) + (h-x)*(y+h)*S(x,y+h))/h**2
1277
1278# Define conversion function for dof's from 2D to scalar label, and
1279# its inverse
1280def ij_to_dof_index(i,j) :
1281 return 33*j+i
1282
1283def inv_ij_to_dof_index(k) :
1284 return [k-33*int(k/33),int(k/33)]
1285
1286
1287# Construct measurement matrix, M, for measurements
1288xs = np.arange(1./14,13./14,1./14); #measurement points
1289
1290M = np.zeros((13,13,33**2));
1291for k in range(33**2) :
1292 c = inv_ij_to_dof_index(k)
1293 for i in range(13) :
1294 for j in range(13) :
1295 M[i,j,k] = phi(xs[i]-h*c[0], xs[j]-h*c[1])
1296M = M.reshape((13**2, 33**2))
1297M = scipy.sparse.csr_matrix(M);
1298
1299# Construct local overlap matrix, A_loc, and identity matrix Id
1300A_loc = np.array([[2./3, -1./6, -1./3, -1./6],
1301 [-1./6, 2./3, -1./6, -1./3],
1302 [-1./3, -1./6, 2./3, -1./6],
1303 [-1./6, -1./3, -1./6, 2./3]])
1304Id = np.eye(33**2,33**2)
1305
1306# Locate boundary labels
1307boundaries = ([ij_to_dof_index(i,0) for i in range(33)] +
1308 [ij_to_dof_index(i,32) for i in range(33)] +
1309 [ij_to_dof_index(0,j+1) for j in range(31)] +
1310 [ij_to_dof_index(32,j+1) for j in range(31)])
1311
1312# Define RHS of FEM linear system, AU = b
1313b = np.ones(33**2)*10*h**2
1314b[boundaries] = 0 #enforce boundary conditions on b
1315
1316
1317
1318
1319
1320###########################################################################
1321###################### forward solver function ############################
1322###########################################################################
1323
1324def forward_solver(theta) :
1325 # Initialize matrix A for FEM linear solve, AU = b
1326 A = np.zeros((33**2,33**2))
1327
1328 # Build A by summing over contribution from each cell
1329 for i in range(32) :
1330 for j in range (32) :
1331 # Find local coefficient in 8x8 grid
1332 theta_loc = theta[int(i/4)+int(j/4)*8]
1333
1334 # Update A by including contribution from cell (i,j)
1335 dof = [ij_to_dof_index(i,j),
1336 ij_to_dof_index(i,j+1),
1337 ij_to_dof_index(i+1,j+1),
1338 ij_to_dof_index(i+1,j)]
1339 A[np.ix_(dof,dof)] += theta_loc * A_loc
1340
1341 # Enforce boundary condition: Zero out rows and columns, then
1342 # put a one back into the diagonal entries.
1343 A[boundaries,:] = 0
1344 A[:,boundaries] = 0
1345 A[boundaries,boundaries] = 1
1346
1347 # Solve linear equation for coefficients, U, and then
1348 # get the Z vector by multiplying by the measurement matrix
1349 u = spsolve(scipy.sparse.csr_matrix(A), b)
1350
1351 z = M * u
1352
1353 return z
1354
1355
1356
1357
1358###########################################################################
1359################# compute log probability, log pi #########################
1360###########################################################################
1361
1362def log_likelihood(theta) :
1363 z = forward_solver(theta)
1364 misfit = z - z_hat
1365 sig = 0.05 #likelihood standard deviation
1366 return -np.dot(misfit,misfit)/(2*sig**2)
1367
1368def log_prior(theta) :
1369 sig_pr = 2 #prior (log) standard deviation
1370 return -np.linalg.norm(np.log(theta))**2/(2*sig_pr**2)
1371
1372def log_posterior(theta) :
1373 return log_likelihood(theta) + log_prior(theta)
1374
1375
1376
1377###########################################################################
1378############# A function to test against known output #####################
1379###########################################################################
1380
1381
1382def verify_against_stored_tests() :
1383 for i in range(10) :
1384 print ("Verifying against data set", i)
1385
1386 # Read the input vector
1387 f_input = open ("../testing/input.{}.txt".format(i), 'r')
1388 theta = np.fromfile(f_input, count=64, sep=" ")
1389
1390 # Then compute both the forward solution and its statistics.
1391 # This is not efficiently written here (it calls the forward
1392 # solver twice), but we don't care about efficiency here since
1393 # we are only computing with ten samples
1394 this_z = forward_solver(theta)
1395 this_log_likelihood = log_likelihood(theta)
1396 this_log_prior = log_prior(theta)
1397
1398 # Then also read the reference output generated by the C++ program:
1399 f_output_z = open ("../testing/output.{}.z.txt".format(i), 'r')
1400 f_output_likelihood = open ("../testing/output.{}.loglikelihood.txt".format(i), 'r')
1401 f_output_prior = open ("../testing/output.{}.logprior.txt".format(i), 'r')
1402
1403 reference_z = np.fromfile(f_output_z, count=13**2, sep=" ")
1404 reference_log_likelihood = float(f_output_likelihood.read())
1405 reference_log_prior = float(f_output_prior.read())
1406
1407 print (" || z-z_ref || : ",
1408 np.linalg.norm(this_z - reference_z))
1409 print (" log likelihood : ",
1410 "Python value=", this_log_likelihood,
1411 "(C++ reference value=", reference_log_likelihood,
1412 ", error=", abs(this_log_likelihood - reference_log_likelihood),
1413 ")")
1414 print (" log prior : ",
1415 "Python value=", this_log_prior,
1416 "(C++ reference value=", reference_log_prior,
1417 ", error=", abs(this_log_prior - reference_log_prior),
1418 ")")
1419
1420
1421def time_forward_solver() :
1422 begin = time.time()
1423
1424 n_runs = 100
1425 for i in range(n_runs) :
1426 # Create a random vector (with entries between 0 and 1), scale
1427 # it by a factor of 4, subtract 2, then take the exponential
1428 # of each entry to get random entries between e^{-2} and
1429 # e^{+2}
1430 theta = np.exp(np.random.rand(64) * 4 - 2)
1431 z = forward_solver(theta)
1432 end = time.time()
1433 print ("Time per forward evaluation:", (end-begin)/n_runs)
1434
1435verify_against_stored_tests()
1436time_forward_solver()
1437@endcode
1438
1439
1440<a name="ann-mcmc-laplace.cc"></a>
1441<h1>Annotated version of mcmc-laplace.cc</h1>
1442 *
1443 *
1444 *
1445 *
1446 * @code
1447 *   /* -----------------------------------------------------------------------------
1448 *   *
1449 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1450 *   * Copyright (C) 2019 by Wolfgang Bangerth
1451 *   *
1452 *   * This file is part of the deal.II code gallery.
1453 *   *
1454 *   * -----------------------------------------------------------------------------
1455 *   *
1456 *   * Author: Wolfgang Bangerth, Colorado State University, 2019.
1457 *   */
1458 *  
1459 *  
1460 *   #include <deal.II/base/timer.h>
1461 *   #include <deal.II/base/multithread_info.h>
1462 *   #include <deal.II/grid/tria.h>
1463 *   #include <deal.II/dofs/dof_handler.h>
1464 *   #include <deal.II/grid/grid_generator.h>
1465 *   #include <deal.II/grid/tria_accessor.h>
1466 *   #include <deal.II/grid/tria_iterator.h>
1467 *   #include <deal.II/dofs/dof_accessor.h>
1468 *   #include <deal.II/fe/fe_q.h>
1469 *   #include <deal.II/dofs/dof_tools.h>
1470 *   #include <deal.II/fe/fe_values.h>
1471 *   #include <deal.II/base/quadrature_lib.h>
1472 *   #include <deal.II/base/function.h>
1473 *   #include <deal.II/numerics/vector_tools.h>
1474 *   #include <deal.II/numerics/matrix_tools.h>
1475 *   #include <deal.II/lac/vector.h>
1476 *   #include <deal.II/lac/full_matrix.h>
1477 *   #include <deal.II/lac/sparse_matrix.h>
1478 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
1479 *   #include <deal.II/lac/solver_cg.h>
1480 *   #include <deal.II/lac/precondition.h>
1481 *   #include <deal.II/lac/sparse_ilu.h>
1482 *  
1483 *   #include <deal.II/numerics/data_out.h>
1484 *  
1485 *   #include <fstream>
1486 *   #include <iostream>
1487 *   #include <random>
1488 *  
1489 *   #include <deal.II/base/logstream.h>
1490 *  
1491 *   using namespace dealii;
1492 *  
1493 *  
1494 * @endcode
1495 *
1496 * The following is a namespace in which we define the solver of the PDE.
1497 * The main class implements an abstract `Interface` class declared at
1498 * the top, which provides for an `evaluate()` function that, given
1499 * a coefficient vector, solves the PDE discussed in the Readme file
1500 * and then evaluates the solution at the 169 mentioned points.
1501 *
1502
1503 *
1504 * The solver follows the basic layout of @ref step_4 "step-4", though it precomputes
1505 * a number of things in the `setup_system()` function, such as the
1506 * evaluation of the matrix that corresponds to the point evaluations,
1507 * as well as the local contributions to matrix and right hand side.
1508 *
1509
1510 *
1511 * Rather than commenting on everything in detail, in the following
1512 * we will only document those things that are not already clear from
1513 * @ref step_4 "step-4" and a small number of other tutorial programs.
1514 *
1515 * @code
1516 *   namespace ForwardSimulator
1517 *   {
1518 *   class Interface
1519 *   {
1520 *   public:
1521 *   virtual Vector<double> evaluate(const Vector<double> &coefficients) = 0;
1522 *  
1523 *   virtual ~Interface() = default;
1524 *   };
1525 *  
1526 *  
1527 *  
1528 *   template <int dim>
1529 *   class PoissonSolver : public Interface
1530 *   {
1531 *   public:
1532 *   PoissonSolver(const unsigned int global_refinements,
1533 *   const unsigned int fe_degree,
1534 *   const std::string &dataset_name);
1535 *   virtual Vector<double>
1536 *   evaluate(const Vector<double> &coefficients) override;
1537 *  
1538 *   private:
1539 *   void make_grid(const unsigned int global_refinements);
1540 *   void setup_system();
1541 *   void assemble_system(const Vector<double> &coefficients);
1542 *   void solve();
1543 *   void output_results(const Vector<double> &coefficients) const;
1544 *  
1545 *   Triangulation<dim> triangulation;
1546 *   FE_Q<dim> fe;
1547 *   DoFHandler<dim> dof_handler;
1548 *  
1549 *   FullMatrix<double> cell_matrix;
1550 *   Vector<double> cell_rhs;
1551 *   std::map<types::global_dof_index,double> boundary_values;
1552 *  
1553 *   SparsityPattern sparsity_pattern;
1554 *   SparseMatrix<double> system_matrix;
1555 *  
1556 *   Vector<double> solution;
1557 *   Vector<double> system_rhs;
1558 *  
1559 *   std::vector<Point<dim>> measurement_points;
1560 *  
1561 *   SparsityPattern measurement_sparsity;
1562 *   SparseMatrix<double> measurement_matrix;
1563 *  
1564 *   TimerOutput timer;
1565 *   unsigned int nth_evaluation;
1566 *  
1567 *   const std::string &dataset_name;
1568 *   };
1569 *  
1570 *  
1571 *  
1572 *   template <int dim>
1573 *   PoissonSolver<dim>::PoissonSolver(const unsigned int global_refinements,
1574 *   const unsigned int fe_degree,
1575 *   const std::string &dataset_name)
1576 *   : fe(fe_degree)
1577 *   , dof_handler(triangulation)
1578 *   , timer(std::cout, TimerOutput::summary, TimerOutput::cpu_times)
1579 *   , nth_evaluation(0)
1580 *   , dataset_name(dataset_name)
1581 *   {
1582 *   make_grid(global_refinements);
1583 *   setup_system();
1584 *   }
1585 *  
1586 *  
1587 *  
1588 *   template <int dim>
1589 *   void PoissonSolver<dim>::make_grid(const unsigned int global_refinements)
1590 *   {
1591 *   Assert(global_refinements >= 3,
1592 *   ExcMessage("This program makes the assumption that the mesh for the "
1593 *   "solution of the PDE is at least as fine as the one used "
1594 *   "in the definition of the coefficient."));
1595 *   GridGenerator::hyper_cube(triangulation, 0, 1);
1596 *   triangulation.refine_global(global_refinements);
1597 *  
1598 *   std::cout << " Number of active cells: " << triangulation.n_active_cells()
1599 *   << std::endl;
1600 *   }
1601 *  
1602 *  
1603 *  
1604 *   template <int dim>
1605 *   void PoissonSolver<dim>::setup_system()
1606 *   {
1607 * @endcode
1608 *
1609 * First define the finite element space:
1610 *
1611 * @code
1612 *   dof_handler.distribute_dofs(fe);
1613 *  
1614 *   std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
1615 *   << std::endl;
1616 *  
1617 * @endcode
1618 *
1619 * Then set up the main data structures that will hold the discrete problem:
1620 *
1621 * @code
1622 *   {
1623 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
1624 *   DoFTools::make_sparsity_pattern(dof_handler, dsp);
1625 *   sparsity_pattern.copy_from(dsp);
1626 *  
1627 *   system_matrix.reinit(sparsity_pattern);
1628 *  
1629 *   solution.reinit(dof_handler.n_dofs());
1630 *   system_rhs.reinit(dof_handler.n_dofs());
1631 *   }
1632 *  
1633 * @endcode
1634 *
1635 * And then define the tools to do point evaluation. We choose
1636 * a set of 13x13 points evenly distributed across the domain:
1637 *
1638 * @code
1639 *   {
1640 *   const unsigned int n_points_per_direction = 13;
1641 *   const double dx = 1. / (n_points_per_direction + 1);
1642 *  
1643 *   for (unsigned int x = 1; x <= n_points_per_direction; ++x)
1644 *   for (unsigned int y = 1; y <= n_points_per_direction; ++y)
1645 *   measurement_points.emplace_back(x * dx, y * dx);
1646 *  
1647 * @endcode
1648 *
1649 * First build a full matrix of the evaluation process. We do this
1650 * even though the matrix is really sparse -- but we don't know
1651 * which entries are nonzero. Later, the `copy_from()` function
1652 * calls build a sparsity pattern and a sparse matrix from
1653 * the dense matrix.
1654 *
1655 * @code
1656 *   Vector<double> weights(dof_handler.n_dofs());
1657 *   FullMatrix<double> full_measurement_matrix(n_points_per_direction *
1658 *   n_points_per_direction,
1659 *   dof_handler.n_dofs());
1660 *  
1661 *   for (unsigned int index = 0; index < measurement_points.size(); ++index)
1662 *   {
1664 *   measurement_points[index],
1665 *   weights);
1666 *   for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
1667 *   full_measurement_matrix(index, i) = weights(i);
1668 *   }
1669 *  
1670 *   measurement_sparsity.copy_from(full_measurement_matrix);
1671 *   measurement_matrix.reinit(measurement_sparsity);
1672 *   measurement_matrix.copy_from(full_measurement_matrix);
1673 *   }
1674 *  
1675 * @endcode
1676 *
1677 * Next build the mapping from cell to the index in the 64-element
1678 * coefficient vector:
1679 *
1680 * @code
1681 *   for (const auto &cell : triangulation.active_cell_iterators())
1682 *   {
1683 *   const unsigned int i = std::floor(cell->center()[0] * 8);
1684 *   const unsigned int j = std::floor(cell->center()[1] * 8);
1685 *  
1686 *   const unsigned int index = i + 8 * j;
1687 *  
1688 *   cell->set_user_index(index);
1689 *   }
1690 *  
1691 * @endcode
1692 *
1693 * Finally prebuild the building blocks of the linear system as
1694 * discussed in the Readme file :
1695 *
1696 * @code
1697 *   {
1698 *   const unsigned int dofs_per_cell = fe.dofs_per_cell;
1699 *  
1700 *   cell_matrix.reinit(dofs_per_cell, dofs_per_cell);
1701 *   cell_rhs.reinit(dofs_per_cell);
1702 *  
1703 *   const QGauss<dim> quadrature_formula(fe.degree+1);
1704 *   const unsigned int n_q_points = quadrature_formula.size();
1705 *  
1706 *   FEValues<dim> fe_values(fe,
1707 *   quadrature_formula,
1709 *   update_JxW_values);
1710 *  
1711 *   fe_values.reinit(dof_handler.begin_active());
1712 *  
1713 *   for (unsigned int q_index = 0; q_index < n_q_points; ++q_index)
1714 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1715 *   {
1716 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1717 *   cell_matrix(i, j) +=
1718 *   (fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
1719 *   fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
1720 *   fe_values.JxW(q_index)); // dx
1721 *  
1722 *   cell_rhs(i) += (fe_values.shape_value(i, q_index) * // phi_i(x_q)
1723 *   10.0 * // f(x_q)
1724 *   fe_values.JxW(q_index)); // dx
1725 *   }
1726 *  
1728 *   0,
1730 *   boundary_values);
1731 *   }
1732 *   }
1733 *  
1734 *  
1735 *  
1736 * @endcode
1737 *
1738 * Given that we have pre-built the matrix and right hand side contributions
1739 * for a (representative) cell, the function that assembles the matrix is
1740 * pretty short and straightforward:
1741 *
1742 * @code
1743 *   template <int dim>
1744 *   void PoissonSolver<dim>::assemble_system(const Vector<double> &coefficients)
1745 *   {
1746 *   Assert(coefficients.size() == 64, ExcInternalError());
1747 *  
1748 *   system_matrix = 0;
1749 *   system_rhs = 0;
1750 *  
1751 *   const unsigned int dofs_per_cell = fe.dofs_per_cell;
1752 *  
1753 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1754 *  
1755 *   for (const auto &cell : dof_handler.active_cell_iterators())
1756 *   {
1757 *   const double coefficient = coefficients(cell->user_index());
1758 *  
1759 *   cell->get_dof_indices(local_dof_indices);
1760 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1761 *   {
1762 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1763 *   system_matrix.add(local_dof_indices[i],
1764 *   local_dof_indices[j],
1765 *   coefficient * cell_matrix(i, j));
1766 *  
1767 *   system_rhs(local_dof_indices[i]) += cell_rhs(i);
1768 *   }
1769 *   }
1770 *  
1771 *   MatrixTools::apply_boundary_values(boundary_values,
1772 *   system_matrix,
1773 *   solution,
1774 *   system_rhs);
1775 *   }
1776 *  
1777 *  
1778 * @endcode
1779 *
1780 * The same is true for the function that solves the linear system:
1781 *
1782 * @code
1783 *   template <int dim>
1784 *   void PoissonSolver<dim>::solve()
1785 *   {
1786 *   SparseILU<double> ilu;
1787 *   ilu.initialize(system_matrix);
1788 *   SolverControl control(100, 1e-10*system_rhs.l2_norm());
1789 *   SolverCG<> solver(control);
1790 *   solver.solve(system_matrix, solution, system_rhs, ilu);
1791 *   }
1792 *  
1793 *  
1794 *  
1795 * @endcode
1796 *
1797 * The following function outputs graphical data for the most recently
1798 * used coefficient and corresponding solution of the PDE. Collecting
1799 * the coefficient values requires translating from the 64-element
1800 * coefficient vector and the cells that correspond to each of these
1801 * elements. The rest remains pretty obvious, with the exception
1802 * of including the number of the current sample into the file name.
1803 *
1804 * @code
1805 *   template <int dim>
1806 *   void
1807 *   PoissonSolver<dim>::output_results(const Vector<double> &coefficients) const
1808 *   {
1809 *   Vector<float> coefficient_values(triangulation.n_active_cells());
1810 *   for (const auto &cell : triangulation.active_cell_iterators())
1811 *   coefficient_values[cell->active_cell_index()] =
1812 *   coefficients(cell->user_index());
1813 *  
1814 *   DataOut<dim> data_out;
1815 *  
1816 *   data_out.attach_dof_handler(dof_handler);
1817 *   data_out.add_data_vector(solution, "solution");
1818 *   data_out.add_data_vector(coefficient_values, "coefficient");
1819 *  
1820 *   data_out.build_patches();
1821 *  
1822 *   std::ofstream output("solution-" +
1823 *   Utilities::int_to_string(nth_evaluation, 10) + ".vtu");
1824 *   data_out.write_vtu(output);
1825 *   }
1826 *  
1827 *  
1828 *  
1829 * @endcode
1830 *
1831 * The following is the main function of this class: Given a coefficient
1832 * vector, it assembles the linear system, solves it, and then evaluates
1833 * the solution at the measurement points by applying the measurement
1834 * matrix to the solution vector. That vector of "measured" values
1835 * is then returned.
1836 *
1837
1838 *
1839 * The function will also output the solution in a graphical format
1840 * if you un-comment the corresponding statement in the third
1841 * code block. However, you may end up with a very large amount
1842 * of data: This code is producing, at the minimum, 10,000 samples
1843 * and creating output for each one of them is surely more data
1844 * than you ever want to see!
1845 *
1846
1847 *
1848 * At the end of the function, we output some timing information
1849 * every 10,000 samples.
1850 *
1851 * @code
1852 *   template <int dim>
1853 *   Vector<double>
1854 *   PoissonSolver<dim>::evaluate(const Vector<double> &coefficients)
1855 *   {
1856 *   {
1857 *   TimerOutput::Scope section(timer, "Building linear systems");
1858 *   assemble_system(coefficients);
1859 *   }
1860 *  
1861 *   {
1862 *   TimerOutput::Scope section(timer, "Solving linear systems");
1863 *   solve();
1864 *   }
1865 *  
1866 *   Vector<double> measurements(measurement_matrix.m());
1867 *   {
1868 *   TimerOutput::Scope section(timer, "Postprocessing");
1869 *  
1870 *   measurement_matrix.vmult(measurements, solution);
1871 *   Assert(measurements.size() == measurement_points.size(),
1872 *   ExcInternalError());
1873 *  
1874 *   /* output_results(coefficients); */
1875 *   }
1876 *  
1877 *   ++nth_evaluation;
1878 *   if (nth_evaluation % 10000 == 0)
1879 *   timer.print_summary();
1880 *  
1881 *   return measurements;
1882 *   }
1883 *   } // namespace ForwardSimulator
1884 *  
1885 *  
1886 * @endcode
1887 *
1888 * The following namespaces define the statistical properties of the Bayesian
1889 * inverse problem. The first is about the definition of the measurement
1890 * statistics (the "likelihood"), which we here assume to be a normal
1891 * distribution @f$N(\mu,\sigma I)@f$ with mean value @f$\mu@f$ given by the
1892 * actual measurement vector (passed as an argument to the constructor
1893 * of the `Gaussian` class and standard deviation @f$\sigma@f$.
1894 *
1895
1896 *
1897 * For reasons of numerical accuracy, it is useful to not return the
1898 * actual likelihood, but its logarithm. This is because these
1899 * values can be very small, occasionally on the order of @f$e^{-100}@f$,
1900 * for which it becomes very difficult to compute accurate
1901 * values.
1902 *
1903 * @code
1904 *   namespace LogLikelihood
1905 *   {
1906 *   class Interface
1907 *   {
1908 *   public:
1909 *   virtual double log_likelihood(const Vector<double> &x) const = 0;
1910 *  
1911 *   virtual ~Interface() = default;
1912 *   };
1913 *  
1914 *  
1915 *   class Gaussian : public Interface
1916 *   {
1917 *   public:
1918 *   Gaussian(const Vector<double> &mu, const double sigma);
1919 *  
1920 *   virtual double log_likelihood(const Vector<double> &x) const override;
1921 *  
1922 *   private:
1923 *   const Vector<double> mu;
1924 *   const double sigma;
1925 *   };
1926 *  
1927 *   Gaussian::Gaussian(const Vector<double> &mu, const double sigma)
1928 *   : mu(mu)
1929 *   , sigma(sigma)
1930 *   {}
1931 *  
1932 *  
1933 *   double Gaussian::log_likelihood(const Vector<double> &x) const
1934 *   {
1935 *   Vector<double> x_minus_mu = x;
1936 *   x_minus_mu -= mu;
1937 *  
1938 *   return -x_minus_mu.norm_sqr() / (2 * sigma * sigma);
1939 *   }
1940 *   } // namespace LogLikelihood
1941 *  
1942 *  
1943 * @endcode
1944 *
1945 * Next up is the "prior" imposed on the coefficients. We assume
1946 * that the logarithms of the entries of the coefficient vector
1947 * are all distributed as a Gaussian with given mean and standard
1948 * deviation. If the logarithms of the coefficients are normally
1949 * distributed, then this implies in particular that the coefficients
1950 * can only be positive, which is a useful property to ensure the
1951 * well-posedness of the forward problem.
1952 *
1953
1954 *
1955 * For the same reasons as for the likelihood above, the interface
1956 * for the prior asks for returning the *logarithm* of the prior,
1957 * instead of the prior probability itself.
1958 *
1959 * @code
1960 *   namespace LogPrior
1961 *   {
1962 *   class Interface
1963 *   {
1964 *   public:
1965 *   virtual double log_prior(const Vector<double> &x) const = 0;
1966 *  
1967 *   virtual ~Interface() = default;
1968 *   };
1969 *  
1970 *  
1971 *   class LogGaussian : public Interface
1972 *   {
1973 *   public:
1974 *   LogGaussian(const double mu, const double sigma);
1975 *  
1976 *   virtual double log_prior(const Vector<double> &x) const override;
1977 *  
1978 *   private:
1979 *   const double mu;
1980 *   const double sigma;
1981 *   };
1982 *  
1983 *   LogGaussian::LogGaussian(const double mu, const double sigma)
1984 *   : mu(mu)
1985 *   , sigma(sigma)
1986 *   {}
1987 *  
1988 *  
1989 *   double LogGaussian::log_prior(const Vector<double> &x) const
1990 *   {
1991 *   double log_of_product = 0;
1992 *  
1993 *   for (const auto &el : x)
1994 *   log_of_product +=
1995 *   -(std::log(el) - mu) * (std::log(el) - mu) / (2 * sigma * sigma);
1996 *  
1997 *   return log_of_product;
1998 *   }
1999 *   } // namespace LogPrior
2000 *  
2001 *  
2002 *  
2003 * @endcode
2004 *
2005 * The Metropolis-Hastings algorithm requires a method to create a new sample
2006 * given a previous sample. We do this by perturbing the current (coefficient)
2007 * sample randomly using a Gaussian distribution centered at the current
2008 * sample. To ensure that the samples' individual entries all remain
2009 * positive, we use a Gaussian distribution in logarithm space -- in other
2010 * words, instead of *adding* a small perturbation with mean value zero,
2011 * we *multiply* the entries of the current sample by a factor that
2012 * is the exponential of a random number with mean zero. (Because the
2013 * exponential of zero is one, this means that the most likely factors
2014 * to multiply the existing sample entries by are close to one. And
2015 * because the exponential of a number is always positive, we never
2016 * get negative samples this way.)
2017 *
2018
2019 *
2020 * But the Metropolis-Hastings sampler doesn't just need a perturbed
2021 * sample @f$y@f$ location given the current sample location @f$x@f$. It also
2022 * needs to know the ratio of the probability of reaching @f$y@f$ from
2023 * @f$x@f$, divided by the probability of reaching @f$x@f$ from @f$y@f$. If we
2024 * were to use a symmetric proposal distribution (e.g., a Gaussian
2025 * distribution centered at @f$x@f$ with a width independent of @f$x@f$), then
2026 * these two probabilities would be the same, and the ratio one. But
2027 * that's not the case for the Gaussian in log space. It's not
2028 * terribly difficult to verify that in that case, for a single
2029 * component the ratio of these probabilities is @f$y_i/x_i@f$, and
2030 * consequently for all components of the vector together, the
2031 * probability is the product of these ratios.
2032 *
2033 * @code
2034 *   namespace ProposalGenerator
2035 *   {
2036 *   class Interface
2037 *   {
2038 *   public:
2039 *   virtual
2040 *   std::pair<Vector<double>,double>
2041 *   perturb(const Vector<double> &current_sample) const = 0;
2042 *  
2043 *   virtual ~Interface() = default;
2044 *   };
2045 *  
2046 *  
2047 *   class LogGaussian : public Interface
2048 *   {
2049 *   public:
2050 *   LogGaussian(const unsigned int random_seed, const double log_sigma);
2051 *  
2052 *   virtual
2053 *   std::pair<Vector<double>,double>
2054 *   perturb(const Vector<double> &current_sample) const override;
2055 *  
2056 *   private:
2057 *   const double log_sigma;
2058 *   mutable std::mt19937 random_number_generator;
2059 *   };
2060 *  
2061 *  
2062 *  
2063 *   LogGaussian::LogGaussian(const unsigned int random_seed,
2064 *   const double log_sigma)
2065 *   : log_sigma(log_sigma)
2066 *   {
2067 *   random_number_generator.seed(random_seed);
2068 *   }
2069 *  
2070 *  
2071 *   std::pair<Vector<double>,double>
2072 *   LogGaussian::perturb(const Vector<double> &current_sample) const
2073 *   {
2074 *   Vector<double> new_sample = current_sample;
2075 *   double product_of_ratios = 1;
2076 *   for (auto &x : new_sample)
2077 *   {
2078 *   const double rnd = std::normal_distribution<>(0, log_sigma)(random_number_generator);
2079 *   const double exp_rnd = std::exp(rnd);
2080 *   x *= exp_rnd;
2081 *   product_of_ratios *= exp_rnd;
2082 *   }
2083 *  
2084 *   return {new_sample, product_of_ratios};
2085 *   }
2086 *  
2087 *   } // namespace ProposalGenerator
2088 *  
2089 *  
2090 * @endcode
2091 *
2092 * The last main class is the Metropolis-Hastings sampler itself.
2093 * If you understand the algorithm behind this method, then
2094 * the following implementation should not be too difficult
2095 * to read. The only thing of relevance is that descriptions
2096 * of the algorithm typically ask whether the *ratio* of two
2097 * probabilities (the "posterior" probabilities of the current
2098 * and the previous samples, where the "posterior" is the product of the
2099 * likelihood and the prior probability) is larger or smaller than a
2100 * randomly drawn number. But because our interfaces return the
2101 * *logarithms* of these probabilities, we now need to take
2102 * the ratio of appropriate exponentials -- which is made numerically
2103 * more stable by considering the exponential of the difference of
2104 * the log probabilities. The only other slight complication is that
2105 * we need to multiply this ratio by the ratio of proposal probabilities
2106 * since we use a non-symmetric proposal distribution.
2107 *
2108
2109 *
2110 * Finally, we note that the output is generated with 7 digits of
2111 * accuracy. (The C++ default is 6 digits.) We do this because,
2112 * as shown in the paper, we can determine the mean value of the
2113 * probability distribution we are sampling here to at least six
2114 * digits of accuracy, and do not want to be limited by the precision
2115 * of the output.
2116 *
2117 * @code
2118 *   namespace Sampler
2119 *   {
2120 *   class MetropolisHastings
2121 *   {
2122 *   public:
2123 *   MetropolisHastings(ForwardSimulator::Interface & simulator,
2124 *   const LogLikelihood::Interface & likelihood,
2125 *   const LogPrior::Interface & prior,
2126 *   const ProposalGenerator::Interface &proposal_generator,
2127 *   const unsigned int random_seed,
2128 *   const std::string & dataset_name);
2129 *  
2130 *   void sample(const Vector<double> &starting_guess,
2131 *   const unsigned int n_samples);
2132 *  
2133 *   private:
2134 *   ForwardSimulator::Interface & simulator;
2135 *   const LogLikelihood::Interface & likelihood;
2136 *   const LogPrior::Interface & prior;
2137 *   const ProposalGenerator::Interface &proposal_generator;
2138 *  
2139 *   std::mt19937 random_number_generator;
2140 *  
2141 *   unsigned int sample_number;
2142 *   unsigned int accepted_sample_number;
2143 *  
2144 *   std::ofstream output_file;
2145 *  
2146 *   void write_sample(const Vector<double> &current_sample,
2147 *   const double current_log_likelihood);
2148 *   };
2149 *  
2150 *  
2151 *   MetropolisHastings::MetropolisHastings(
2152 *   ForwardSimulator::Interface & simulator,
2153 *   const LogLikelihood::Interface & likelihood,
2154 *   const LogPrior::Interface & prior,
2155 *   const ProposalGenerator::Interface &proposal_generator,
2156 *   const unsigned int random_seed,
2157 *   const std::string & dataset_name)
2158 *   : simulator(simulator)
2159 *   , likelihood(likelihood)
2160 *   , prior(prior)
2161 *   , proposal_generator(proposal_generator)
2162 *   , sample_number(0)
2163 *   , accepted_sample_number(0)
2164 *   {
2165 *   output_file.open("samples-" + dataset_name + ".txt");
2166 *   output_file.precision(7);
2167 *  
2168 *   random_number_generator.seed(random_seed);
2169 *   }
2170 *  
2171 *  
2172 *   void MetropolisHastings::sample(const Vector<double> &starting_guess,
2173 *   const unsigned int n_samples)
2174 *   {
2175 *   std::uniform_real_distribution<> uniform_distribution(0, 1);
2176 *  
2177 *   Vector<double> current_sample = starting_guess;
2178 *   double current_log_posterior =
2179 *   (likelihood.log_likelihood(simulator.evaluate(current_sample)) +
2180 *   prior.log_prior(current_sample));
2181 *  
2182 *   ++sample_number;
2183 *   ++accepted_sample_number;
2184 *   write_sample(current_sample, current_log_posterior);
2185 *  
2186 *   for (unsigned int k = 1; k < n_samples; ++k, ++sample_number)
2187 *   {
2188 *   std::pair<Vector<double>,double>
2189 *   perturbation = proposal_generator.perturb(current_sample);
2190 *   const Vector<double> trial_sample = std::move (perturbation.first);
2191 *   const double perturbation_probability_ratio = perturbation.second;
2192 *  
2193 *   const double trial_log_posterior =
2194 *   (likelihood.log_likelihood(simulator.evaluate(trial_sample)) +
2195 *   prior.log_prior(trial_sample));
2196 *  
2197 *   if (std::exp(trial_log_posterior - current_log_posterior) * perturbation_probability_ratio
2198 *   >=
2199 *   uniform_distribution(random_number_generator))
2200 *   {
2201 *   current_sample = trial_sample;
2202 *   current_log_posterior = trial_log_posterior;
2203 *  
2204 *   ++accepted_sample_number;
2205 *   }
2206 *  
2207 *   write_sample(current_sample, current_log_posterior);
2208 *   }
2209 *   }
2210 *  
2211 *  
2212 *  
2213 *   void MetropolisHastings::write_sample(const Vector<double> &current_sample,
2214 *   const double current_log_posterior)
2215 *   {
2216 *   output_file << current_log_posterior << '\t';
2217 *   output_file << accepted_sample_number << '\t';
2218 *   for (const auto &x : current_sample)
2219 *   output_file << x << ' ';
2220 *   output_file << '\n';
2221 *  
2222 *   output_file.flush();
2223 *   }
2224 *   } // namespace Sampler
2225 *  
2226 *  
2227 * @endcode
2228 *
2229 * The final function is `main()`, which simply puts all of these pieces
2230 * together into one. The "exact solution", i.e., the "measurement values"
2231 * we use for this program are tabulated to make it easier for other
2232 * people to use in their own implementations of this benchmark. These
2233 * values created using the same main class above, but using 8 mesh
2234 * refinements and using a Q3 element -- i.e., using a much more accurate
2235 * method than the one we use in the forward simulator for generating
2236 * samples below (which uses 5 global mesh refinement steps and a Q1
2237 * element). If you wanted to regenerate this set of numbers, then
2238 * the following code snippet would do that:
2239 * <div class=CodeFragmentInTutorialComment>
2240 * @code
2241 * /* Set the exact coefficient: */
2242 * Vector<double> exact_coefficients(64);
2243 * for (auto &el : exact_coefficients)
2244 * el = 1.;
2245 * exact_coefficients(9) = exact_coefficients(10) = exact_coefficients(17) =
2246 * exact_coefficients(18) = 0.1;
2247 * exact_coefficients(45) = exact_coefficients(46) = exact_coefficients(53) =
2248 * exact_coefficients(54) = 10.;
2249 *
2250
2251 * /* Compute the "correct" solution vector: */
2252 * const Vector<double> exact_solution =
2253 * ForwardSimulator::PoissonSolver<2>(/* global_refinements = */ 8,
2254 * /* fe_degree = */ 3,
2255 * /* prefix = */ "exact")
2256 * .evaluate(exact_coefficients);
2257 * @endcode
2258 * </div>
2259 *
2260 * @code
2261 *   int main()
2262 *   {
2263 *   const bool testing = true;
2264 *  
2265 * @endcode
2266 *
2267 * Run with one thread, so as to not step on other processes
2268 * doing the same at the same time. It turns out that the problem
2269 * is also so small that running with more than one thread
2270 * *increases* the runtime.
2271 *
2272 * @code
2274 *  
2275 *   const unsigned int random_seed = (testing ? 1U : std::random_device()());
2276 *   const std::string dataset_name = Utilities::to_string(random_seed, 10);
2277 *  
2278 *   const Vector<double> exact_solution(
2279 *   { 0.06076511762259369, 0.09601910120848481,
2280 *   0.1238852517838584, 0.1495184117375201,
2281 *   0.1841596127549784, 0.2174525028261122,
2282 *   0.2250996160898698, 0.2197954769002993,
2283 *   0.2074695698370926, 0.1889996477663016,
2284 *   0.1632722532153726, 0.1276782480038186,
2285 *   0.07711845915789312, 0.09601910120848552,
2286 *   0.2000589533367983, 0.3385592591951766,
2287 *   0.3934300024647806, 0.4040223892461541,
2288 *   0.4122329537843092, 0.4100480091545554,
2289 *   0.3949151637189968, 0.3697873264791232,
2290 *   0.33401826235924, 0.2850397806663382,
2291 *   0.2184260032478671, 0.1271121156350957,
2292 *   0.1238852517838611, 0.3385592591951819,
2293 *   0.7119285162766475, 0.8175712861756428,
2294 *   0.6836254116578105, 0.5779452419831157,
2295 *   0.5555615956136897, 0.5285181561736719,
2296 *   0.491439702849224, 0.4409367494853282,
2297 *   0.3730060082060772, 0.2821694983395214,
2298 *   0.1610176733857739, 0.1495184117375257,
2299 *   0.3934300024647929, 0.8175712861756562,
2300 *   0.9439154625527653, 0.8015904115095128,
2301 *   0.6859683749254024, 0.6561235366960599,
2302 *   0.6213197201867315, 0.5753611315000049,
2303 *   0.5140091754526823, 0.4325325506354165,
2304 *   0.3248315148915482, 0.1834600412730086,
2305 *   0.1841596127549917, 0.4040223892461832,
2306 *   0.6836254116578439, 0.8015904115095396,
2307 *   0.7870119561144977, 0.7373108331395808,
2308 *   0.7116558878070463, 0.6745179049094283,
2309 *   0.6235300574156917, 0.5559332704045935,
2310 *   0.4670304994474178, 0.3499809143811,
2311 *   0.19688263746294, 0.2174525028261253,
2312 *   0.4122329537843404, 0.5779452419831566,
2313 *   0.6859683749254372, 0.7373108331396063,
2314 *   0.7458811983178246, 0.7278968022406559,
2315 *   0.6904793535357751, 0.6369176452710288,
2316 *   0.5677443693743215, 0.4784738764865867,
2317 *   0.3602190632823262, 0.2031792054737325,
2318 *   0.2250996160898818, 0.4100480091545787,
2319 *   0.5555615956137137, 0.6561235366960938,
2320 *   0.7116558878070715, 0.727896802240657,
2321 *   0.7121928678670187, 0.6712187391428729,
2322 *   0.6139157775591492, 0.5478251665295381,
2323 *   0.4677122687599031, 0.3587654911000848,
2324 *   0.2050734291675918, 0.2197954769003094,
2325 *   0.3949151637190157, 0.5285181561736911,
2326 *   0.6213197201867471, 0.6745179049094407,
2327 *   0.690479353535786, 0.6712187391428787,
2328 *   0.6178408289359514, 0.5453605027237883,
2329 *   0.489575966490909, 0.4341716881061278,
2330 *   0.3534389974779456, 0.2083227496961347,
2331 *   0.207469569837099, 0.3697873264791366,
2332 *   0.4914397028492412, 0.5753611315000203,
2333 *   0.6235300574157017, 0.6369176452710497,
2334 *   0.6139157775591579, 0.5453605027237935,
2335 *   0.4336604929612851, 0.4109641743019312,
2336 *   0.3881864790111245, 0.3642640090182592,
2337 *   0.2179599909280145, 0.1889996477663011,
2338 *   0.3340182623592461, 0.4409367494853381,
2339 *   0.5140091754526943, 0.5559332704045969,
2340 *   0.5677443693743304, 0.5478251665295453,
2341 *   0.4895759664908982, 0.4109641743019171,
2342 *   0.395727260284338, 0.3778949322004734,
2343 *   0.3596268271857124, 0.2191250268948948,
2344 *   0.1632722532153683, 0.2850397806663325,
2345 *   0.373006008206081, 0.4325325506354207,
2346 *   0.4670304994474315, 0.4784738764866023,
2347 *   0.4677122687599041, 0.4341716881061055,
2348 *   0.388186479011099, 0.3778949322004602,
2349 *   0.3633362567187364, 0.3464457261905399,
2350 *   0.2096362321365655, 0.1276782480038148,
2351 *   0.2184260032478634, 0.2821694983395252,
2352 *   0.3248315148915535, 0.3499809143811097,
2353 *   0.3602190632823333, 0.3587654911000799,
2354 *   0.3534389974779268, 0.3642640090182283,
2355 *   0.35962682718569, 0.3464457261905295,
2356 *   0.3260728953424643, 0.180670595355394,
2357 *   0.07711845915789244, 0.1271121156350963,
2358 *   0.1610176733857757, 0.1834600412730144,
2359 *   0.1968826374629443, 0.2031792054737354,
2360 *   0.2050734291675885, 0.2083227496961245,
2361 *   0.2179599909279998, 0.2191250268948822,
2362 *   0.2096362321365551, 0.1806705953553887,
2363 *   0.1067965550010013 });
2364 *  
2365 * @endcode
2366 *
2367 * Now run the forward simulator for samples:
2368 *
2369 * @code
2370 *   ForwardSimulator::PoissonSolver<2> laplace_problem(
2371 *   /* global_refinements = */ 5,
2372 *   /* fe_degree = */ 1,
2373 *   dataset_name);
2374 *   LogLikelihood::Gaussian log_likelihood(exact_solution, 0.05);
2375 *   LogPrior::LogGaussian log_prior(0, 2);
2376 *   ProposalGenerator::LogGaussian proposal_generator(
2377 *   random_seed, 0.09); /* so that the acceptance ratio is ~0.24 */
2378 *   Sampler::MetropolisHastings sampler(laplace_problem,
2379 *   log_likelihood,
2380 *   log_prior,
2381 *   proposal_generator,
2382 *   random_seed,
2383 *   dataset_name);
2384 *  
2385 *   Vector<double> starting_coefficients(64);
2386 *   for (auto &el : starting_coefficients)
2387 *   el = 1.;
2388 *   sampler.sample(starting_coefficients,
2389 *   (testing ? 250 * 40 /* takes 10 seconds */
2390 *   :
2391 *   100000000 /* takes 1.5 days */
2392 *   ));
2393 *   }
2394 * @endcode
2395
2396
2397*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
static void set_thread_limit(const unsigned int max_threads=numbers::invalid_unsigned_int)
void initialize(const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
real_type norm_sqr() const
Point< 2 > first
Definition grid_out.cc:4623
#define Assert(cond, exc)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
std::vector< index_type > data
Definition mpi.cc:735
std::size_t size
Definition mpi.cc:734
Expression floor(const Expression &x)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ diagonal
Matrix is diagonal.
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:479
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
void create_point_source_vector(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const Point< spacedim, double > &p, Vector< double > &rhs_vector)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
STL namespace.
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation