DOpE
statpdeproblem.h
Go to the documentation of this file.
1 
24 #ifndef STAT_PDE_PROBLEM_H_
25 #define STAT_PDE_PROBLEM_H_
26 
27 #include "pdeprobleminterface.h"
28 #include "integrator.h"
29 #include "parameterreader.h"
30 #include "statevector.h"
31 #include "stateproblem.h"
32 
33 #include <lac/vector.h>
34 
35 #include <lac/block_sparsity_pattern.h>
36 #include <lac/block_sparse_matrix.h>
37 
38 #include "pdeproblemcontainer.h"
39 #include "pdeinterface.h"
40 #include "functionalinterface.h"
41 #include "dirichletdatainterface.h"
42 #include "dopeexception.h"
43 #include "newtonsolver.h"
44 #include "cglinearsolver.h"
45 #include "gmreslinearsolver.h"
46 #include "directlinearsolver.h"
47 #include "solutionextractor.h"
48 
49 #include <base/data_out_base.h>
50 #include <base/utilities.h>
51 #include <numerics/data_out.h>
52 #include <numerics/matrix_tools.h>
53 #include <numerics/vector_tools.h>
54 #include <base/function.h>
55 #include <lac/sparse_matrix.h>
56 #include <lac/compressed_simple_sparsity_pattern.h>
57 #include <lac/block_sparsity_pattern.h>
58 #include <lac/sparse_direct.h>
59 
60 #include <fstream>
61 namespace DOpE
62 {
77  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
78  typename VECTOR, int dealdim>
79  class StatPDEProblem : public PDEProblemInterface<PROBLEM, VECTOR, dealdim>
80  {
81  public:
94  template<typename INTEGRATORDATACONT>
95  StatPDEProblem(PROBLEM *OP, DOpEtypes::VectorStorageType state_behavior,
96  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
97  int base_priority = 0);
98 
111  template<typename INTEGRATORDATACONT>
112  StatPDEProblem(PROBLEM *OP, DOpEtypes::VectorStorageType state_behavior,
113  ParameterReader &param_reader, INTEGRATORDATACONT& idc1,
114  INTEGRATORDATACONT& idc2, int base_priority = 0);
115 
116  virtual
117  ~StatPDEProblem();
118 
119  /******************************************************/
120 
126  static void
127  declare_params(ParameterReader &param_reader);
128 
129  /******************************************************/
130 
136  void
137  ReInit();
138 
139  /******************************************************/
140 
146  void
148 
149  /******************************************************/
150 
156  void
158 
159  /******************************************************/
160 
178  template<class DWRC, class PDE>
179  void
180  ComputeRefinementIndicators(DWRC& dwrc, PDE& pde);
181 
182  /******************************************************/
183 
189  void
190  StateSizeInfo(std::stringstream& out)
191  {
192  GetU().PrintInfos(out);
193  }
194 
195  /******************************************************/
196 
207  void
208  WriteToFile(const VECTOR &v, std::string name, std::string outfile,
209  std::string dof_type, std::string filetype);
210 
211  /******************************************************/
212 
226  void
227  WriteToFileElementwise(const Vector<double> &v, std::string name,
228  std::string outfile, std::string dof_type, std::string filetype);
229 
230  /******************************************************/
231 
240  void
241  WriteToFile(const ControlVector<VECTOR> &v, std::string name, std::string dof_type);
242 
250  void
251  WriteToFile(const std::vector<double> &/*v*/,
252  std::string /*outfile*/)
253  {
254  abort();
255  }
256 
257  protected:
261  void
263 
264  /******************************************************/
265 
275  void
277 
278  /******************************************************/
282  const StateVector<VECTOR> &
283  GetU() const
284  {
285  return u_;
286  }
289  {
290  return u_;
291  }
295  const StateVector<VECTOR> &
296  GetZForEE() const
297  {
298  return z_for_ee_;
299  }
302  {
303  return z_for_ee_;
304  }
305 
306  NONLINEARSOLVER&
307  GetNonlinearSolver(std::string type);
308  INTEGRATOR&
310  {
311  return integrator_;
312  }
313 
314  const bool&
316  {
317  return build_state_matrix_;
318  }
319 
320  const bool&
322  {
323  return build_adjoint_matrix_;
324  }
325 
326  const bool&
328  {
329  return state_reinit_;
330  }
331 
332  const bool&
334  {
335  return adjoint_reinit_;
336  }
337 
338 
339  bool&
341  {
342  return build_state_matrix_;
343  }
344 
345  bool&
347  {
348  return build_adjoint_matrix_;
349  }
350 
351  bool&
353  {
354  return state_reinit_;
355  }
356 
357  bool&
359  {
360  return adjoint_reinit_;
361  }
362 
363 
364  private:
369  void
370  AddUDD()
371  {
372  for (auto it = this->GetUserDomainData().begin();
373  it != this->GetUserDomainData().end(); it++)
374  {
375  this->GetIntegrator().AddDomainData(it->first, it->second);
376  }
377  }
378 
383  void
384  DeleteUDD()
385  {
386  for (auto it = this->GetUserDomainData().begin();
387  it != this->GetUserDomainData().end(); it++)
388  {
389  this->GetIntegrator().DeleteDomainData(it->first);
390  }
391  }
392 
393  StateVector<VECTOR> u_;
394  StateVector<VECTOR> z_for_ee_;
395 
396  INTEGRATOR integrator_;
397  NONLINEARSOLVER nonlinear_state_solver_;
398  NONLINEARSOLVER nonlinear_adjoint_solver_;
399 
400  bool build_state_matrix_;
401  bool build_adjoint_matrix_;
402  bool state_reinit_;
403  bool adjoint_reinit_;
404 
405  int n_patches_;
406 
407  friend class SolutionExtractor<
408  StatPDEProblem<NONLINEARSOLVER, INTEGRATOR, PROBLEM, VECTOR, dealdim>,
409  VECTOR> ;
410  };
411 
412  /*************************************************************************/
413  /*****************************IMPLEMENTATION******************************/
414  /*************************************************************************/
415  using namespace dealii;
416 
417  /******************************************************/
418  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
419  typename VECTOR, int dealdim>
420  void
422  ParameterReader &param_reader)
423  {
424  NONLINEARSOLVER::declare_params(param_reader);
425  param_reader.SetSubsection("output parameters");
426  param_reader.declare_entry("number of patches", "0",
427  Patterns::Integer(0));
428 
429  }
430  /******************************************************/
431 
432  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
433  typename VECTOR, int dealdim>
434  template<typename INTEGRATORDATACONT>
436  PROBLEM *OP, DOpEtypes::VectorStorageType state_behavior,
437  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
438  int base_priority)
439  : PDEProblemInterface<PROBLEM, VECTOR, dealdim>(OP, base_priority), u_(
440  OP->GetSpaceTimeHandler(), state_behavior, param_reader), z_for_ee_(
441  OP->GetSpaceTimeHandler(), state_behavior, param_reader), integrator_(
442  idc), nonlinear_state_solver_(integrator_, param_reader), nonlinear_adjoint_solver_(
443  integrator_, param_reader)
444  {
445  param_reader.SetSubsection("output parameters");
446  n_patches_ = param_reader.get_integer("number of patches");
447  //PDEProblems should be ReInited
448  {
449  state_reinit_ = true;
450  adjoint_reinit_ = true;
451  }
452  }
453 
454  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
455  typename VECTOR, int dealdim>
456  template<typename INTEGRATORDATACONT>
458  PROBLEM *OP, DOpEtypes::VectorStorageType state_behavior,
459  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
460  INTEGRATORDATACONT& idc2, int base_priority)
461  : PDEProblemInterface<PROBLEM, VECTOR, dealdim>(OP, base_priority), u_(
462  OP->GetSpaceTimeHandler(), state_behavior, param_reader), z_for_ee_(
463  OP->GetSpaceTimeHandler(), state_behavior, param_reader), integrator_(
464  idc, idc2), nonlinear_state_solver_(integrator_, param_reader), nonlinear_adjoint_solver_(
465  integrator_, param_reader)
466  {
467  //PDEProblems should be ReInited
468  {
469  state_reinit_ = true;
470  adjoint_reinit_ = true;
471  }
472  }
473 
474  /******************************************************/
475 
476  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
477  typename VECTOR, int dealdim>
479  {
480  }
481 
482  /******************************************************/
483 
484  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
485  typename VECTOR, int dealdim>
486  NONLINEARSOLVER&
488  std::string type)
489  {
490  if (type == "state")
491  {
492  return nonlinear_state_solver_;
493  }
494  else if (type == "adjoint_for_ee")
495  {
496  return nonlinear_adjoint_solver_;
497  }
498  else
499  {
500  throw DOpEException("No Solver for Problem type:`" + type + "' found",
501  "StatReducedProblem::GetNonlinearSolver");
502 
503  }
504  }
505 
506  /******************************************************/
507 
508  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
509  typename VECTOR, int dealdim>
510  void
512  {
514 
515  //Some Solvers must be reinited when called
516  // Better we have subproblems, so that solver can be reinited here
517  {
518  state_reinit_ = true;
519  adjoint_reinit_ = true;
520  }
521 
522  build_state_matrix_ = true;
523  build_adjoint_matrix_ = true;
524 
525  GetU().ReInit();
526  GetZForEE().ReInit();
527  }
528 
529  /******************************************************/
530 
531  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
532  typename VECTOR, int dealdim>
533  void
535  {
536  this->InitializeFunctionalValues(this->GetProblem()->GetNFunctionals());
537  this->GetOutputHandler()->Write("Computing State Solution:",
538  4 + this->GetBasePriority());
539 
540  this->SetProblemType("state");
541  auto& problem = this->GetProblem()->GetStateProblem();
542  if (state_reinit_ == true)
543  {
544  GetNonlinearSolver("state").ReInit(problem);
545  state_reinit_ = false;
546  }
547 
548  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
549  AddUDD();
550  build_state_matrix_ = this->GetNonlinearSolver("state").NonlinearSolve(
551  problem, (GetU().GetSpacialVector()), true, build_state_matrix_);
552  DeleteUDD();
553 
554  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
555 
556  this->GetOutputHandler()->Write((GetU().GetSpacialVector()),
557  "State" + this->GetPostIndex(), problem.GetDoFType());
558 
559  }
560 
561  /******************************************************/
562 
563  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
564  typename VECTOR, int dealdim>
565  void
567  DOpEtypes::WeightComputation weight_comp)
568  {
569  this->GetOutputHandler()->Write("Computing Dual for Error Estimation:",
570  4 + this->GetBasePriority());
571  if (weight_comp == DOpEtypes::higher_order_interpolation)
572  {
573  this->SetProblemType("adjoint_for_ee");
574  }
575  else
576  {
577  throw DOpEException("Unknown WeightComputation",
578  "StatPDEProblem::ComputeDualForErrorEstimation");
579  }
580 
581  // auto& problem = this->GetProblem()->GetStateProblem();//Hier ist adjoint problem einzufuegen
582  auto& problem = *(this->GetProblem());
583  if (adjoint_reinit_ == true)
584  {
585  GetNonlinearSolver("adjoint_for_ee").ReInit(problem);
586  adjoint_reinit_ = false;
587  }
588 
589  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
590 
591  this->GetIntegrator().AddDomainData("state",
592  &(GetU().GetSpacialVector()));
593  AddUDD();
594 
595  build_adjoint_matrix_ =
596  this->GetNonlinearSolver("adjoint_for_ee").NonlinearSolve(problem,
597  (GetZForEE().GetSpacialVector()), true, build_adjoint_matrix_);
598 
599  this->GetIntegrator().DeleteDomainData("state");
600  DeleteUDD();
601 
602  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
603 
604  this->GetOutputHandler()->Write((GetZForEE().GetSpacialVector()),
605  "Adjoint_for_ee" + this->GetPostIndex(), problem.GetDoFType());
606 
607  }
608 
609  /******************************************************/
610 
611  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
612  typename VECTOR, int dealdim>
613  void
615  {
616  this->ComputeReducedState();
617 
618  this->GetOutputHandler()->Write("Computing Functionals:",
619  4 + this->GetBasePriority());
620 
621  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
622 
623  this->GetIntegrator().AddDomainData("state",
624  &(GetU().GetSpacialVector()));
625  AddUDD();
626 
627  for (unsigned int i = 0; i < this->GetProblem()->GetNFunctionals(); i++)
628  {
629  double ret = 0;
630  bool found = false;
631 
632  this->SetProblemType("aux_functional", i);
633  if (this->GetProblem()->GetFunctionalType().find("domain")
634  != std::string::npos)
635  {
636  found = true;
637  ret += this->GetIntegrator().ComputeDomainScalar(
638  *(this->GetProblem()));
639  }
640  if (this->GetProblem()->GetFunctionalType().find("point")
641  != std::string::npos)
642  {
643  found = true;
644  ret += this->GetIntegrator().ComputePointScalar(
645  *(this->GetProblem()));
646  }
647  if (this->GetProblem()->GetFunctionalType().find("boundary")
648  != std::string::npos)
649  {
650  found = true;
651  ret += this->GetIntegrator().ComputeBoundaryScalar(
652  *(this->GetProblem()));
653  }
654  if (this->GetProblem()->GetFunctionalType().find("face")
655  != std::string::npos)
656  {
657  found = true;
658  ret += this->GetIntegrator().ComputeFaceScalar(*(this->GetProblem()));
659  }
660 
661  if (!found)
662  {
663  throw DOpEException(
664  "Unknown Functional Type: "
665  + this->GetProblem()->GetFunctionalType(),
666  "StatPDEProblem::ComputeReducedFunctionals");
667  }
668  this->GetFunctionalValues()[i].push_back(ret);
669  std::stringstream out;
670  this->GetOutputHandler()->InitOut(out);
671  out << this->GetProblem()->GetFunctionalName() << ": " << ret;
672  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
673  }
674 
675  this->GetIntegrator().DeleteDomainData("state");
676  DeleteUDD();
677 
678  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
679  }
680 
681  /******************************************************/
682 
683  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
684  typename VECTOR, int dealdim>
685  void
687  const StateVector<VECTOR>& st_u)
688  {
689  this->InitializeFunctionalValues(this->GetProblem()->GetNFunctionals());
690  GetU() = st_u;
691 
692  this->GetOutputHandler()->Write("Computing Functionals:",
693  4 + this->GetBasePriority());
694 
695  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
696 
697  this->GetIntegrator().AddDomainData("state",
698  &(GetU().GetSpacialVector()));
699  AddUDD();
700 
701  for (unsigned int i = 0; i < this->GetProblem()->GetNFunctionals(); i++)
702  {
703  double ret = 0;
704  bool found = false;
705 
706  this->SetProblemType("aux_functional", i);
707  if (this->GetProblem()->GetFunctionalType().find("domain")
708  != std::string::npos)
709  {
710  found = true;
711  ret += this->GetIntegrator().ComputeDomainScalar(
712  *(this->GetProblem()));
713  }
714  if (this->GetProblem()->GetFunctionalType().find("point")
715  != std::string::npos)
716  {
717  found = true;
718  ret += this->GetIntegrator().ComputePointScalar(
719  *(this->GetProblem()));
720  }
721  if (this->GetProblem()->GetFunctionalType().find("boundary")
722  != std::string::npos)
723  {
724  found = true;
725  ret += this->GetIntegrator().ComputeBoundaryScalar(
726  *(this->GetProblem()));
727  }
728  if (this->GetProblem()->GetFunctionalType().find("face")
729  != std::string::npos)
730  {
731  found = true;
732  ret += this->GetIntegrator().ComputeFaceScalar(*(this->GetProblem()));
733  }
734 
735  if (!found)
736  {
737  throw DOpEException(
738  "Unknown Functional Type: "
739  + this->GetProblem()->GetFunctionalType(),
740  "StatPDEProblem::ComputeReducedFunctionals");
741  }
742  this->GetFunctionalValues()[i].push_back(ret);
743  std::stringstream out;
744  this->GetOutputHandler()->InitOut(out);
745  out << this->GetProblem()->GetFunctionalName() << ": " << ret;
746  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
747  }
748 
749  this->GetIntegrator().DeleteDomainData("state");
750  DeleteUDD();
751 
752  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
753  }
754 
755  /******************************************************/
756 
757  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
758  typename VECTOR, int dealdim>
759  template<class DWRC, class PDE>
760  void
762  DWRC& dwrc, PDE& pde)
763  {
764  //Attach the ResidualModifier to the PDE.
765  pde.ResidualModifier = boost::bind<void>(
766  boost::mem_fn(&DWRC::ResidualModifier), boost::ref(dwrc), _1);
767  pde.VectorResidualModifier = boost::bind<void>(
768  boost::mem_fn(&DWRC::VectorResidualModifier), boost::ref(dwrc), _1);
769  //first we reinit the dwrdatacontainer (this
770  //sets the weight-vectors to their correct length)
771  const unsigned int n_elements =
772  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler().get_tria().n_active_cells();
773  dwrc.ReInit(n_elements);
774  //If we need the dual solution, compute it
775  if (dwrc.NeedDual())
776  this->ComputeDualForErrorEstimation(dwrc.GetWeightComputation());
777 
778  //some output
779  this->GetOutputHandler()->Write("Computing Error Indicators:",
780  4 + this->GetBasePriority());
781 
782  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
783 
784  //add the primal and (if needed) dual solution to the integrator
785  this->GetIntegrator().AddDomainData("state",
786  &(GetU().GetSpacialVector()));
787  if (dwrc.NeedDual())
788  this->GetIntegrator().AddDomainData("adjoint_for_ee",
789  &(GetZForEE().GetSpacialVector()));
790  AddUDD();
791 
792  this->SetProblemType("error_evaluation");
793 
794  //prepare the weights...
795  dwrc.PrepareWeights(GetU(), GetZForEE());
796 
797  //now we finally compute the refinement indicators
798  this->GetIntegrator().ComputeRefinementIndicators(*this->GetProblem(),
799  dwrc);
800 
801  // release the lock on the refinement indicators (see dwrcontainer.h)
802  dwrc.ReleaseLock();
803 
804  const float error = dwrc.GetError();
805 
806  // clear the data
807  dwrc.ClearWeightData();
808  this->GetIntegrator().DeleteDomainData("state");
809  if (dwrc.NeedDual())
810  this->GetIntegrator().DeleteDomainData("adjoint_for_ee");
811  DeleteUDD();
812  this->GetProblem()->DeleteAuxiliaryFromIntegrator(
813  this->GetIntegrator());
814 
815  std::stringstream out;
816  this->GetOutputHandler()->InitOut(out);
817  out << "Error estimate using " << dwrc.GetName();
818  if (dwrc.NeedDual())
819  out << " for the " << this->GetProblem()->GetFunctionalName();
820  out << ": " << error;
821  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
822  this->GetOutputHandler()->WriteElementwise(dwrc.GetErrorIndicators(),
823  "Error_Indicators" + this->GetPostIndex(),
824  this->GetProblem()->GetDoFType());
825 
826  }
827 
828  /******************************************************/
829 
830  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
831  typename VECTOR, int dealdim>
832  void
834  const Vector<double> &v, std::string name, std::string outfile,
835  std::string dof_type, std::string filetype)
836  {
837  if (dof_type == "state")
838  {
839  auto& data_out =
840  this->GetProblem()->GetSpaceTimeHandler()->GetDataOut();
841  data_out.attach_dof_handler(
842  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler());
843 
844  data_out.add_data_vector(v, name);
845  data_out.build_patches(n_patches_);
846 
847  std::ofstream output(outfile.c_str());
848 
849  if (filetype == ".vtk")
850  {
851  data_out.write_vtk(output);
852  }
853  else
854  {
855  throw DOpEException(
856  "Don't know how to write filetype `" + filetype + "'!",
857  "StatPDEProblem::WriteToFileElementwise");
858  }
859  data_out.clear();
860  }
861  else
862  {
863  throw DOpEException("No such DoFHandler `" + dof_type + "'!",
864  "StatPDEProblem::WriteToFileElementwise");
865  }
866  }
867 
868  /******************************************************/
869 
870  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
871  typename VECTOR, int dealdim>
872  void
874  const VECTOR &v, std::string name, std::string outfile,
875  std::string dof_type, std::string filetype)
876  {
877  if (dof_type == "state")
878  {
879  auto& data_out =
880  this->GetProblem()->GetSpaceTimeHandler()->GetDataOut();
881  data_out.attach_dof_handler(
882  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler());
883 
884  data_out.add_data_vector(v, name);
885  //TODO: mapping[0] is a workaround, as deal does not support interpolate
886  // boundary_values with a mapping collection at this point.
887  data_out.build_patches(
888  this->GetProblem()->GetSpaceTimeHandler()->GetMapping()[0]);
889 
890  std::ofstream output(outfile.c_str());
891 
892  if (filetype == ".vtk")
893  {
894  data_out.write_vtk(output);
895  }
896  else if (filetype == ".gpl")
897  {
898  data_out.write_gnuplot(output);
899  }
900  else
901  {
902  throw DOpEException(
903  "Don't know how to write filetype `" + filetype + "'!",
904  "StatPDEProblem::WriteToFile");
905  }
906  data_out.clear();
907  }
908  else
909  {
910  throw DOpEException("No such DoFHandler `" + dof_type + "'!",
911  "StatPDEProblem::WriteToFile");
912  }
913  }
914 
915  /******************************************************/
916 
917  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
918  typename VECTOR, int dealdim>
919  void
921  const ControlVector<VECTOR> &v, std::string name, std::string dof_type)
922  {
923  throw DOpEException("This Problem does not support ControlVectors","StatPDEProblem::WriteToFile");
924  }
925 }
926 #endif
StateVector< VECTOR > & GetU()
Definition: statpdeproblem.h:288
void StateSizeInfo(std::stringstream &out)
Definition: statpdeproblem.h:190
void WriteToFileElementwise(const Vector< double > &v, std::string name, std::string outfile, std::string dof_type, std::string filetype)
Definition: statpdeproblem.h:833
Definition: statpdeproblem.h:79
void ComputeReducedFunctionals()
Definition: statpdeproblem.h:614
void WriteToFile(const std::vector< double > &, std::string)
Definition: statpdeproblem.h:251
Definition: parameterreader.h:36
const bool & GetBuildAdjointMatrix() const
Definition: statpdeproblem.h:321
INTEGRATOR & GetIntegrator()
Definition: statpdeproblem.h:309
const StateVector< VECTOR > & GetZForEE() const
Definition: statpdeproblem.h:296
const bool & GetStateReinit() const
Definition: statpdeproblem.h:327
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Definition: parameterreader.h:98
const bool & GetBuildStateMatrix() const
Definition: statpdeproblem.h:315
NONLINEARSOLVER & GetNonlinearSolver(std::string type)
Definition: statpdeproblem.h:487
Definition: controlvector.h:49
Definition: pdeprobleminterface.h:53
StateVector< VECTOR > & GetZForEE()
Definition: statpdeproblem.h:301
bool & GetBuildAdjointMatrix()
Definition: statpdeproblem.h:346
static void declare_params(ParameterReader &param_reader)
Definition: statpdeproblem.h:421
Definition: statevector.h:50
bool & GetStateReinit()
Definition: statpdeproblem.h:352
const StateVector< VECTOR > & GetU() const
Definition: statpdeproblem.h:283
bool & GetAdjointReinit()
Definition: statpdeproblem.h:358
void WriteToFile(const VECTOR &v, std::string name, std::string outfile, std::string dof_type, std::string filetype)
Definition: statpdeproblem.h:873
virtual ~StatPDEProblem()
Definition: statpdeproblem.h:478
void ComputeRefinementIndicators(DWRC &dwrc, PDE &pde)
Definition: statpdeproblem.h:761
int get_integer(const std::string &entry_name)
Definition: parameterreader.h:126
VectorStorageType
Definition: dopetypes.h:120
void ComputeDualForErrorEstimation(DOpEtypes::WeightComputation)
Definition: statpdeproblem.h:566
void SetSubsection(const std::string subsection)
Definition: parameterreader.h:93
void ReInit()
Definition: statpdeproblem.h:511
WeightComputation
Definition: dopetypes.h:80
const std::map< std::string, const VECTOR * > & GetUserDomainData() const
Definition: reducedprobleminterface.h:301
bool & GetBuildStateMatrix()
Definition: statpdeproblem.h:340
StatPDEProblem(PROBLEM *OP, DOpEtypes::VectorStorageType state_behavior, ParameterReader &param_reader, INTEGRATORDATACONT &idc, int base_priority=0)
Definition: statpdeproblem.h:435
const bool & GetAdjointReinit() const
Definition: statpdeproblem.h:333
void ComputeReducedState()
Definition: statpdeproblem.h:534
Definition: dopeexception.h:35