[Trilinos-Users] 2 moocho questions.

Bartlett, Roscoe A rabartl at sandia.gov
Thu Jul 24 10:24:49 MDT 2008


________________________________
From: trilinos-users-bounces at software.sandia.gov [mailto:trilinos-users-bounces at software.sandia.gov] On Behalf Of Luca Heltai
Sent: Thursday, July 24, 2008 7:52 AM
To: trilinos-users at software.sandia.gov
Subject: [Trilinos-Users] 2 moocho questions.

Dear Trilinos users,

I have a couple of questions about moocho....

I would like to keep track of what is happening during the minimization process by making some "realtime" plots. I'm using the EpetraExt Model Evaluator interface, and I'm letting the package compute finite differences for the objective DgDp.

First question: is there a way to find out from within the EvalModel call the iteration number of the moocho solver?

I'd like to place some code as

static int last_iteration = -1;

if(last_iteration != get_current_iteration_here() ) {
    plot_data();
    last_iteration = get_current_iteration_here();
}


in order to have the intermediate solutions plotted at each iteration step. (each iteration in my problem is about a minute long, so this will not impact performances too much...).
The MOOCHO base class IterationPack::AlgorithmTracker<file:///C:/_mystuff/PROJECTS/Trilinos.base/Trilinos/packages/moocho/browser/doc/html/classIterationPack_1_1AlgorithmTracker.html> is designed allow the kind of observations that you want.  See the subclass in the file:

    Trilinos/packages/moocho/src/MoochoPack/src/std/MoochoPack_MoochoTrackerConsoleStd.[hpp,cpp]

for how to get at all of this data.  This "tracker" object gets called at the end of every optimization iteration and you have clean access to all stored algorithm data.

You should copy the files MoochoPack_MoochoTrackerConsoleStd.[hpp,cpp], change the name of the subclass, and make it do what you want.

You can then set your subclass tracker object with something like:
MoochoThyraSolver moochoThyraSolver;
moochoThyraSolver.getSolver().set_track(
  Teuchos::rcp(new LucasGreatAlgorithmTracker(...))
   );
I don't know if I have any concrete examples that show this but you could check (i.e. grep on 'set_track' and see what comes up).

Let me know if you have any problems with this.

Second question: is there a way to pass the finite difference calculator an Epetra_CrsGraph so that only the existing entries are calculated? The system I'm solving is pretty big and blockwise sparse, so this would really help a lot...
The way things are set up right now, you have to provide d(f)/d(x) yourself as an Epetra_CrsMatrix and the finite-difference capability computes d(f)/d(p) for you as well as O(np) directional derivatives with g(x,p) to get the reduced gradient d(g_hat)/d(p).

Are you asking for the sparse finite-difference computation of d(f)/d(p) as a sparse Epetra_CrsMatrix object?  If you have a large design space (i.e. np is larger), then you will want to use an adjoint-based optimization method.  To get this, you will have to provide fully derivatives yourself for d(f)/d(p), d(g)/d(x), and d(g)/d(p).

Long term, you should look into using automatic differentiation with Sacado to compute the derivatives you need.  We don't have any concrete examples in Trilinos set up to do this yet but we will at some point.

These and more details can be found at:

    http://trilinos.sandia.gov/packages/docs/r8.0/packages/moocho/doc/html/index.html

Cheers,

- Ross
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20080724/3bf7ecc2/attachment.html 


More information about the Trilinos-Users mailing list