00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include "LOCA_MultiPredictor_Tangent.H"
00043 #include "LOCA_GlobalData.H"
00044 #include "NOX_Utils.H"
00045 #include "LOCA_ErrorCheck.H"
00046 #include "LOCA_MultiContinuation_AbstractGroup.H"
00047 #include "LOCA_MultiContinuation_ExtendedGroup.H"
00048 #include "LOCA_MultiContinuation_ExtendedVector.H"
00049 #include "LOCA_MultiContinuation_ExtendedMultiVector.H"
00050 #include "Teuchos_ParameterList.hpp"
00051
00052 LOCA::MultiPredictor::Tangent::Tangent(
00053 const Teuchos::RCP<LOCA::GlobalData>& global_data,
00054 const Teuchos::RCP<Teuchos::ParameterList>& predParams,
00055 const Teuchos::RCP<Teuchos::ParameterList>& solverParams) :
00056 globalData(global_data),
00057 linSolverParams(solverParams),
00058 fdfdp(),
00059 tangent(),
00060 secant(),
00061 initialized(false)
00062 {
00063 }
00064
00065 LOCA::MultiPredictor::Tangent::~Tangent()
00066 {
00067 }
00068
00069 LOCA::MultiPredictor::Tangent::Tangent(
00070 const LOCA::MultiPredictor::Tangent& source,
00071 NOX::CopyType type) :
00072 globalData(source.globalData),
00073 linSolverParams(source.linSolverParams),
00074 fdfdp(),
00075 tangent(),
00076 secant(),
00077 initialized(source.initialized)
00078 {
00079 if (source.initialized) {
00080 fdfdp = source.fdfdp->clone(type);
00081
00082 tangent = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedMultiVector>(source.tangent->clone(type));
00083
00084 secant = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedVector>(source.secant->clone(type));
00085 }
00086 }
00087
00088 LOCA::MultiPredictor::AbstractStrategy&
00089 LOCA::MultiPredictor::Tangent::operator=(
00090 const LOCA::MultiPredictor::AbstractStrategy& s)
00091 {
00092 const LOCA::MultiPredictor::Tangent& source =
00093 dynamic_cast<const LOCA::MultiPredictor::Tangent&>(s);
00094
00095 if (this != &source) {
00096 globalData = source.globalData;
00097 linSolverParams = source.linSolverParams;
00098 initialized = source.initialized;
00099
00100 if (source.initialized) {
00101 fdfdp = source.fdfdp->clone(NOX::DeepCopy);
00102
00103 tangent = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedMultiVector>(source.tangent->clone(NOX::DeepCopy));
00104
00105 secant = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedVector>(source.secant->clone(NOX::DeepCopy));
00106 }
00107 }
00108
00109 return *this;
00110 }
00111
00112 Teuchos::RCP<LOCA::MultiPredictor::AbstractStrategy>
00113 LOCA::MultiPredictor::Tangent::clone(NOX::CopyType type) const
00114 {
00115 return Teuchos::rcp(new Tangent(*this, type));
00116 }
00117
00118 NOX::Abstract::Group::ReturnType
00119 LOCA::MultiPredictor::Tangent::compute(
00120 bool baseOnSecant, const vector<double>& stepSize,
00121 LOCA::MultiContinuation::ExtendedGroup& grp,
00122 const LOCA::MultiContinuation::ExtendedVector& prevXVec,
00123 const LOCA::MultiContinuation::ExtendedVector& xVec)
00124 {
00125 string callingFunction = "LOCA::MultiPredictor::Tangent::compute()";
00126 NOX::Abstract::Group::ReturnType status, finalStatus;
00127
00128 if (globalData->locaUtils->isPrintType(NOX::Utils::StepperDetails))
00129 globalData->locaUtils->out() <<
00130 "\n\tCalling Predictor with method: Tangent" << std::endl;
00131
00132
00133 int numParams = stepSize.size();
00134
00135
00136 Teuchos::RCP<LOCA::MultiContinuation::AbstractGroup> underlyingGroup
00137 = grp.getUnderlyingGroup();
00138
00139 if (!initialized) {
00140
00141
00142 fdfdp = underlyingGroup->getX().createMultiVector(numParams+1,
00143 NOX::ShapeCopy);
00144
00145
00146 tangent = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedMultiVector>(xVec.createMultiVector(numParams, NOX::ShapeCopy));
00147
00148
00149 secant = Teuchos::rcp_dynamic_cast<LOCA::MultiContinuation::ExtendedVector>(xVec.clone(NOX::ShapeCopy));
00150
00151 initialized = true;
00152 }
00153
00154
00155 Teuchos::RCP<NOX::Abstract::MultiVector> tanX =
00156 tangent->getXMultiVec();
00157 Teuchos::RCP<NOX::Abstract::MultiVector::DenseMatrix> tanP =
00158 tangent->getScalars();
00159
00160
00161 const vector<int>& conParamIDs = grp.getContinuationParameterIDs();
00162
00163
00164 finalStatus = underlyingGroup->computeDfDpMulti(conParamIDs, *fdfdp, false);
00165 globalData->locaErrorCheck->checkReturnType(finalStatus, callingFunction);
00166
00167 vector<int> index_dfdp(conParamIDs.size());
00168 for (unsigned int i=0; i<conParamIDs.size(); i++)
00169 index_dfdp[i] = i+1;
00170 Teuchos::RCP<NOX::Abstract::MultiVector>dfdp =
00171 fdfdp->subView(index_dfdp);
00172
00173
00174 for (unsigned int i=0; i<conParamIDs.size(); i++)
00175 (*dfdp)[i].scale(-1.0);
00176
00177
00178 status = underlyingGroup->computeJacobian();
00179 finalStatus =
00180 globalData->locaErrorCheck->combineAndCheckReturnTypes(status, finalStatus,
00181 callingFunction);
00182
00183
00184 status = underlyingGroup->applyJacobianInverseMultiVector(*linSolverParams,
00185 *dfdp,
00186 *tanX);
00187 finalStatus =
00188 globalData->locaErrorCheck->combineAndCheckReturnTypes(status, finalStatus,
00189 callingFunction);
00190
00191
00192 tanP->putScalar(0.0);
00193 for (unsigned int i=0; i<conParamIDs.size(); i++)
00194 (*tanP)(i,i) = 1.0;
00195
00196
00197 setPredictorOrientation(baseOnSecant, stepSize, grp, prevXVec,
00198 xVec, *secant, *tangent);
00199
00200 return finalStatus;
00201 }
00202
00203 NOX::Abstract::Group::ReturnType
00204 LOCA::MultiPredictor::Tangent::evaluate(
00205 const vector<double>& stepSize,
00206 const LOCA::MultiContinuation::ExtendedVector& xVec,
00207 LOCA::MultiContinuation::ExtendedMultiVector& result) const
00208 {
00209
00210 int numParams = stepSize.size();
00211
00212 for (int i=0; i<numParams; i++)
00213 result[i].update(1.0, xVec, stepSize[i], (*tangent)[i], 0.0);
00214
00215 return NOX::Abstract::Group::Ok;
00216 }
00217
00218 NOX::Abstract::Group::ReturnType
00219 LOCA::MultiPredictor::Tangent::computeTangent(
00220 LOCA::MultiContinuation::ExtendedMultiVector& v)
00221 {
00222 v = *tangent;
00223
00224 return NOX::Abstract::Group::Ok;
00225 }
00226
00227 bool
00228 LOCA::MultiPredictor::Tangent::isTangentScalable() const
00229 {
00230 return true;
00231 }