00001 #include "EpetraExt_MultiPointModelEvaluator.h"
00002
00003 EpetraExt::MultiPointModelEvaluator::MultiPointModelEvaluator(
00004 Teuchos::RefCountPtr<EpetraExt::ModelEvaluator> underlyingME_,
00005 const Teuchos::RefCountPtr<EpetraExt::MultiMpiComm> &globalComm_,
00006 const std::vector<Epetra_Vector*> initGuessVec_,
00007 Teuchos::RefCountPtr<std::vector< Teuchos::RefCountPtr<Epetra_Vector> > > q_vec_,
00008 Teuchos::RefCountPtr<std::vector< Teuchos::RefCountPtr<Epetra_Vector> > > matching_vec_
00009 ) :
00010 underlyingME(underlyingME_),
00011 globalComm(globalComm_),
00012 timeStepsOnTimeDomain(globalComm_->NumTimeStepsOnDomain()),
00013 numTimeDomains(globalComm_->NumSubDomains()),
00014 timeDomain(globalComm_->SubDomainRank()),
00015 rowStencil(0),
00016 rowIndex(0),
00017 underlyingNg(0),
00018 q_vec(q_vec_),
00019 matching_vec(matching_vec_)
00020 {
00021 if (globalComm->MyPID()==0) {
00022 cout << "----------MultiPoint Partition Info------------"
00023 << "\n\tNumProcs = " << globalComm->NumProc()
00024 << "\n\tSpatial Decomposition = " << globalComm->SubDomainComm().NumProc()
00025 << "\n\tNumber of Domains = " << numTimeDomains
00026 << "\n\tSteps on Domain 0 = " << timeStepsOnTimeDomain
00027 << "\n\tTotal Number of Steps = " << globalComm->NumTimeSteps();
00028 cout << "\n-----------------------------------------------" << endl;
00029 }
00030
00031
00032
00033
00034 split_W = Teuchos::rcp_dynamic_cast<Epetra_RowMatrix>(underlyingME->create_W());
00035
00036 rowStencil = new std::vector< std::vector<int> >(timeStepsOnTimeDomain);
00037 rowIndex = new std::vector<int>;
00038 for (int i=0; i < timeStepsOnTimeDomain; i++) {
00039 (*rowStencil)[i].push_back(0);
00040 (*rowIndex).push_back(i + globalComm->FirstTimeStepOnDomain());
00041 }
00042
00043 block_W = Teuchos::rcp(new EpetraExt::BlockCrsMatrix(*split_W,
00044 *rowStencil, *rowIndex, *globalComm));
00045
00046
00047 EpetraExt::ModelEvaluator::OutArgs underlyingOutArgs = underlyingME->createOutArgs();
00048
00049 underlyingNg = underlyingOutArgs.Ng();
00050 if (underlyingNg) {
00051 if (underlyingOutArgs.supports(OUT_ARG_DgDp,0,0).supports(DERIV_TRANS_MV_BY_ROW))
00052 orientation_DgDp = DERIV_TRANS_MV_BY_ROW;
00053 else
00054 orientation_DgDp = DERIV_MV_BY_COL;
00055 }
00056
00057
00058 TEST_FOR_EXCEPT(underlyingOutArgs.Np()!=2);
00059
00060
00061 const Epetra_Map& split_map = split_W->RowMatrixRowMap();
00062 num_p0 = underlyingME_->get_p_map(0)->NumMyElements();
00063 if (underlyingNg) num_g0 = underlyingME_->get_g_map(0)->NumMyElements();
00064 else num_g0 = 0;
00065 num_dg0dp0 = num_g0 * num_p0;
00066
00067
00068 block_x = new EpetraExt::BlockVector(split_map, block_W->RowMap());
00069 block_f = new EpetraExt::BlockVector(*block_x);
00070 block_DfDp = new EpetraExt::BlockMultiVector(split_map, block_W->RowMap(), num_p0);
00071 if (underlyingNg)
00072 block_DgDx = new EpetraExt::BlockMultiVector(split_map, block_W->RowMap(), num_g0);
00073
00074
00075 split_x = Teuchos::rcp(new Epetra_Vector(split_map));
00076 split_f = Teuchos::rcp(new Epetra_Vector(split_map));
00077 split_DfDp = Teuchos::rcp(new Epetra_MultiVector(split_map, num_p0));
00078 if (underlyingNg)
00079 split_DgDx = Teuchos::rcp(new Epetra_MultiVector(split_map, num_g0));
00080 if (underlyingNg) {
00081 if(orientation_DgDp == DERIV_TRANS_MV_BY_ROW)
00082 split_DgDp = Teuchos::rcp(new Epetra_MultiVector(*(underlyingME_->get_p_map(0)), num_g0));
00083 else
00084 split_DgDp = Teuchos::rcp(new Epetra_MultiVector(*(underlyingME_->get_g_map(0)), num_p0));
00085 }
00086 if (underlyingNg)
00087 split_g = Teuchos::rcp(new Epetra_Vector(*(underlyingME_->get_g_map(0))));
00088
00089
00090 derivMV_DfDp = new EpetraExt::ModelEvaluator::DerivativeMultiVector(split_DfDp);
00091 deriv_DfDp = new EpetraExt::ModelEvaluator::Derivative(*derivMV_DfDp);
00092 if (underlyingNg) {
00093 derivMV_DgDx = new EpetraExt::ModelEvaluator::DerivativeMultiVector(split_DgDx, DERIV_TRANS_MV_BY_ROW);
00094 deriv_DgDx = new EpetraExt::ModelEvaluator::Derivative(*derivMV_DgDx);
00095 derivMV_DgDp = new EpetraExt::ModelEvaluator::DerivativeMultiVector(split_DgDp, orientation_DgDp);
00096 deriv_DgDp = new EpetraExt::ModelEvaluator::Derivative(*derivMV_DgDp);
00097 }
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107 solution_init = Teuchos::rcp(new EpetraExt::BlockVector(*block_x));
00108 for (int i=0; i < timeStepsOnTimeDomain; i++)
00109 solution_init->LoadBlockValues(*(initGuessVec_[i]), (*rowIndex)[i]);
00110
00111
00112
00113 if (Teuchos::is_null(matching_vec)) matchingProblem = false;
00114 else matchingProblem = true;
00115
00116 if (matchingProblem) {
00117 TEST_FOR_EXCEPT(matching_vec->size()!=timeStepsOnTimeDomain);
00118 TEST_FOR_EXCEPT(!(*matching_vec)[0]->Map().SameAs(*(underlyingME_->get_g_map(0))));
00119 TEST_FOR_EXCEPT(num_g0 != 1);
00120 }
00121 }
00122
00123 EpetraExt::MultiPointModelEvaluator::~MultiPointModelEvaluator()
00124 {
00125 delete block_x;
00126 delete block_f;
00127 delete block_DfDp;
00128 if (underlyingNg) delete block_DgDx;
00129 delete rowStencil;
00130 delete rowIndex;
00131
00132 delete derivMV_DfDp;
00133 delete deriv_DfDp;
00134 if (underlyingNg) {
00135 delete derivMV_DgDx;
00136 delete deriv_DgDx;
00137 delete derivMV_DgDp;
00138 delete deriv_DgDp;
00139 }
00140 }
00141
00142 Teuchos::RefCountPtr<const Epetra_Map> EpetraExt::MultiPointModelEvaluator::get_x_map() const
00143 {
00144 return Teuchos::rcp(&(block_W->OperatorDomainMap()), false);
00145 }
00146
00147 Teuchos::RefCountPtr<const Epetra_Map> EpetraExt::MultiPointModelEvaluator::get_f_map() const
00148 {
00149 return get_x_map();
00150 }
00151
00152 Teuchos::RefCountPtr<const Epetra_Map> EpetraExt::MultiPointModelEvaluator::get_p_map(int l) const
00153 {
00154 return underlyingME->get_p_map(l);
00155 }
00156
00157 Teuchos::RefCountPtr<const Epetra_Map> EpetraExt::MultiPointModelEvaluator::get_g_map(int j) const
00158 {
00159 return underlyingME->get_g_map(j);
00160 }
00161
00162 Teuchos::RefCountPtr<const Epetra_Vector> EpetraExt::MultiPointModelEvaluator::get_x_init() const
00163 {
00164 return solution_init;
00165 }
00166
00167 Teuchos::RefCountPtr<const Epetra_Vector> EpetraExt::MultiPointModelEvaluator::get_p_init(int l) const
00168 {
00169 return underlyingME->get_p_init(l);
00170 }
00171
00172 Teuchos::RefCountPtr<Epetra_Operator> EpetraExt::MultiPointModelEvaluator::create_W() const
00173 {
00174 return block_W;
00175 }
00176
00177 EpetraExt::ModelEvaluator::InArgs EpetraExt::MultiPointModelEvaluator::createInArgs() const
00178 {
00179
00180 InArgsSetup inArgs;
00181 inArgs.setModelEvalDescription(this->description());
00182 inArgs.set_Np(1);
00183 inArgs.setSupports(IN_ARG_x,true);
00184 return inArgs;
00185 }
00186
00187 EpetraExt::ModelEvaluator::OutArgs EpetraExt::MultiPointModelEvaluator::createOutArgs() const
00188 {
00189
00190 OutArgsSetup outArgs;
00191 outArgs.setModelEvalDescription(this->description());
00192 outArgs.set_Np_Ng(1, underlyingNg);
00193 outArgs.setSupports(OUT_ARG_f,true);
00194 outArgs.setSupports(OUT_ARG_W,true);
00195 outArgs.set_W_properties(
00196 DerivativeProperties(
00197 DERIV_LINEARITY_NONCONST
00198 ,DERIV_RANK_FULL
00199 ,true
00200 )
00201 );
00202 outArgs.setSupports(OUT_ARG_DfDp,0,DERIV_MV_BY_COL);
00203 outArgs.set_DfDp_properties(
00204 0,DerivativeProperties(
00205 DERIV_LINEARITY_CONST
00206 ,DERIV_RANK_DEFICIENT
00207 ,true
00208 )
00209 );
00210
00211 if (underlyingNg) {
00212 outArgs.setSupports(OUT_ARG_DgDx,0,DERIV_TRANS_MV_BY_ROW);
00213 outArgs.set_DgDx_properties(
00214 0,DerivativeProperties(
00215 DERIV_LINEARITY_NONCONST
00216 ,DERIV_RANK_DEFICIENT
00217 ,true
00218 )
00219 );
00220 outArgs.setSupports(OUT_ARG_DgDp,0,0, orientation_DgDp);
00221 outArgs.set_DgDp_properties(
00222 0,0,DerivativeProperties(
00223 DERIV_LINEARITY_NONCONST
00224 ,DERIV_RANK_DEFICIENT
00225 ,true
00226 )
00227 );
00228 }
00229 return outArgs;
00230 }
00231
00232 void EpetraExt::MultiPointModelEvaluator::evalModel( const InArgs& inArgs,
00233 const OutArgs& outArgs ) const
00234 {
00235
00236 EpetraExt::ModelEvaluator::InArgs underlyingInArgs = underlyingME->createInArgs();
00237 EpetraExt::ModelEvaluator::OutArgs underlyingOutArgs = underlyingME->createOutArgs();
00238
00239
00240
00241
00242
00243
00244
00245
00246 Teuchos::RefCountPtr<const Epetra_Vector> p_in = inArgs.get_p(0);
00247 if (p_in.get()) underlyingInArgs.set_p(0, p_in);
00248
00249 Teuchos::RefCountPtr<const Epetra_Vector> x_in = inArgs.get_x();
00250 block_x->Epetra_Vector::operator=(*x_in);
00251
00252
00253 Teuchos::RefCountPtr<Epetra_Vector> f_out = outArgs.get_f();
00254
00255 Teuchos::RefCountPtr<Epetra_Operator> W_out = outArgs.get_W();
00256 Teuchos::RefCountPtr<EpetraExt::BlockCrsMatrix> W_block =
00257 Teuchos::rcp_dynamic_cast<EpetraExt::BlockCrsMatrix>(W_out);
00258
00259 Teuchos::RefCountPtr<Epetra_Vector> g_out;
00260 if (underlyingNg) g_out = outArgs.get_g(0);
00261 if (g_out.get()) g_out->PutScalar(0.0);
00262
00263 EpetraExt::ModelEvaluator::Derivative DfDp_out = outArgs.get_DfDp(0);
00264
00265 EpetraExt::ModelEvaluator::Derivative DgDx_out;
00266 EpetraExt::ModelEvaluator::Derivative DgDp_out;
00267 if (underlyingNg) {
00268 DgDx_out = outArgs.get_DgDx(0);
00269 DgDp_out = outArgs.get_DgDp(0,0);
00270 if (!DgDx_out.isEmpty()) DgDx_out.getMultiVector()->PutScalar(0.0);
00271 if (!DgDp_out.isEmpty()) DgDp_out.getMultiVector()->PutScalar(0.0);
00272 }
00273
00274
00275
00276 bool need_g = g_out.get();
00277 if (matchingProblem)
00278 if ( !DgDx_out.isEmpty() || !DgDp_out.isEmpty() ) need_g = true;
00279
00280
00281
00282 for (int i=0; i < timeStepsOnTimeDomain; i++) {
00283
00284
00285 underlyingInArgs.set_p(1, (*q_vec)[i]);
00286
00287
00288 block_x->ExtractBlockValues(*split_x, (*rowIndex)[i]);
00289 underlyingInArgs.set_x(split_x);
00290
00291
00292 if (f_out.get()) underlyingOutArgs.set_f(split_f);
00293
00294 if (need_g) underlyingOutArgs.set_g(0, split_g);
00295
00296 if (W_out.get()) underlyingOutArgs.set_W(split_W);
00297
00298 if (!DfDp_out.isEmpty()) underlyingOutArgs.set_DfDp(0, *deriv_DfDp);
00299
00300 if (!DgDx_out.isEmpty()) underlyingOutArgs.set_DgDx(0, *deriv_DgDx);
00301
00302 if (!DgDp_out.isEmpty()) underlyingOutArgs.set_DgDp(0, 0, *deriv_DgDp);
00303
00304
00305 underlyingME->evalModel(underlyingInArgs, underlyingOutArgs);
00306
00307
00308
00309 if (matchingProblem) {
00310 if (need_g) {
00311 double diff = (*split_g)[0] - (*(*matching_vec)[i])[0];
00312 double nrmlz = fabs((*(*matching_vec)[i])[0]) + 1.0e-6;
00313 (*split_g)[0] = 0.5 * diff * diff/(nrmlz*nrmlz);
00314 if (!DgDx_out.isEmpty()) split_DgDx->Scale(diff/(nrmlz*nrmlz));
00315 if (!DgDp_out.isEmpty()) split_DgDp->Scale(diff/(nrmlz*nrmlz));
00316 }
00317 }
00318
00319
00320 if (f_out.get()) block_f->LoadBlockValues(*split_f, (*rowIndex)[i]);
00321 if (W_out.get()) W_block->LoadBlock(*split_W, i, 0);
00322
00323 if (!DfDp_out.isEmpty()) block_DfDp->LoadBlockValues(*split_DfDp, (*rowIndex)[i]);
00324 if (!DgDx_out.isEmpty()) block_DgDx->LoadBlockValues(*split_DgDx, (*rowIndex)[i]);
00325
00326
00327 if (g_out.get()) g_out->Update(1.0, *split_g, 1.0);
00328
00329 if (!DgDp_out.isEmpty())
00330 DgDp_out.getMultiVector()->Update(1.0, *split_DgDp, 1.0);
00331
00332 }
00333
00334
00335 if (f_out.get()) f_out->operator=(*block_f);
00336 if (!DfDp_out.isEmpty())
00337 DfDp_out.getMultiVector()->operator=(*block_DfDp);
00338 if (!DgDx_out.isEmpty())
00339 DgDx_out.getMultiVector()->operator=(*block_DgDx);
00340
00341
00342 if (numTimeDomains > 1) {
00343 double factorToZeroOutCopies = 0.0;
00344 if (globalComm->SubDomainComm().MyPID()==0) factorToZeroOutCopies = 1.0;
00345 if (g_out.get()) {
00346 (*g_out).Scale(factorToZeroOutCopies);
00347 double* vPtr = &((*g_out)[0]);
00348 Epetra_Vector tmp = *(g_out.get());
00349 globalComm->SumAll( &(tmp[0]), vPtr, num_g0);
00350 }
00351 if (!DgDp_out.isEmpty()) {
00352 DgDp_out.getMultiVector()->Scale(factorToZeroOutCopies);
00353 double* mvPtr = (*DgDp_out.getMultiVector())[0];
00354 Epetra_MultiVector tmp = *(DgDp_out.getMultiVector());
00355 globalComm->SumAll(tmp[0], mvPtr, num_dg0dp0);
00356 }
00357 }
00358 }