Go to the documentation of this file.00001
00002
00003
00004
00005
00006 #include "MOR_ReducedOrderModelFactory.hpp"
00007
00008 #include "MOR_ReducedSpaceFactory.hpp"
00009 #include "MOR_ReducedSpace.hpp"
00010
00011 #include "MOR_BasisOps.hpp"
00012
00013 #include "MOR_SampleDofListFactory.hpp"
00014
00015 #include "MOR_ReducedOrderModelEvaluator.hpp"
00016 #include "MOR_PetrovGalerkinOperatorFactory.hpp"
00017 #include "MOR_GaussNewtonOperatorFactory.hpp"
00018
00019 #include "MOR_ContainerUtils.hpp"
00020
00021 #include "Teuchos_Tuple.hpp"
00022 #include "Teuchos_TestForException.hpp"
00023
00024 #include <string>
00025 #include <stdexcept>
00026
00027 namespace MOR {
00028
00029 using ::Teuchos::RCP;
00030 using ::Teuchos::rcp;
00031 using ::Teuchos::nonnull;
00032 using ::Teuchos::ParameterList;
00033 using ::Teuchos::sublist;
00034 using ::Teuchos::Tuple;
00035 using ::Teuchos::tuple;
00036
00037 ReducedOrderModelFactory::ReducedOrderModelFactory(
00038 const Teuchos::RCP<ReducedSpaceFactory> &spaceFactory,
00039 const RCP<ParameterList> &parentParams) :
00040 spaceFactory_(spaceFactory),
00041 params_(extractModelOrderReductionParams(parentParams))
00042 {
00043
00044 }
00045
00046 RCP<ParameterList> ReducedOrderModelFactory::extractModelOrderReductionParams(const RCP<ParameterList> ¶ms)
00047 {
00048 return sublist(params, "Model Order Reduction");
00049 }
00050
00051 RCP<ParameterList> ReducedOrderModelFactory::extractReducedOrderModelParams(const RCP<ParameterList> ¶ms)
00052 {
00053 return sublist(params, "Reduced-Order Model");
00054 }
00055
00056 RCP<EpetraExt::ModelEvaluator> ReducedOrderModelFactory::create(const RCP<EpetraExt::ModelEvaluator> &child)
00057 {
00058 RCP<EpetraExt::ModelEvaluator> result = child;
00059
00060 if (useReducedOrderModel()) {
00061 const RCP<ParameterList> romParams = extractReducedOrderModelParams(params_);
00062
00063 const Tuple<std::string, 2> allowedProjectionTypes = tuple<std::string>("Galerkin Projection", "Minimum Residual");
00064 const std::string projectionType = romParams->get("System Reduction", allowedProjectionTypes[0]);
00065 TEUCHOS_TEST_FOR_EXCEPTION(!contains(allowedProjectionTypes, projectionType),
00066 std::out_of_range,
00067 projectionType + " not in " + allowedProjectionTypes.toString());
00068
00069 const RCP<const ReducedSpace> reducedSpace = spaceFactory_->create(romParams);
00070 const RCP<const Epetra_MultiVector> basis = spaceFactory_->getBasis(romParams);
00071
00072 if (projectionType == allowedProjectionTypes[0]) {
00073 const RCP<const Epetra_MultiVector> projector = spaceFactory_->getProjector(romParams);
00074 const RCP<ReducedOperatorFactory> opFactory(new PetrovGalerkinOperatorFactory(basis, projector));
00075 result = rcp(new ReducedOrderModelEvaluator(child, reducedSpace, opFactory));
00076 } else if (projectionType == allowedProjectionTypes[1]) {
00077 RCP<ReducedOperatorFactory> opFactory;
00078
00079 const RCP<const Epetra_Operator> collocationOperator =
00080 spaceFactory_->getSamplingOperator(romParams, *child->get_x_map());
00081 if (nonnull(collocationOperator)) {
00082 opFactory = rcp(new GaussNewtonMetricOperatorFactory(basis, collocationOperator));
00083 } else {
00084 opFactory = rcp(new GaussNewtonOperatorFactory(basis));
00085 }
00086
00087 result = rcp(new ReducedOrderModelEvaluator(child, reducedSpace, opFactory));
00088 } else {
00089 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Should not happen");
00090 }
00091 }
00092
00093 return result;
00094 }
00095
00096 bool ReducedOrderModelFactory::useReducedOrderModel() const
00097 {
00098 return extractReducedOrderModelParams(params_)->get("Activate", false);
00099 }
00100
00101 }