1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
/* *****************************************************************
    MESQUITE -- The Mesh Quality Improvement Toolkit

    Copyright 2004 Sandia Corporation and Argonne National
    Laboratory.  Under the terms of Contract DE-AC04-94AL85000
    with Sandia Corporation, the U.S. Government retains certain
    rights in this software.

    This library is free software; you can redistribute it and/or
    modify it under the terms of the GNU Lesser General Public
    License as published by the Free Software Foundation; either
    version 2.1 of the License, or (at your option) any later version.

    This library is distributed in the hope that it will be useful,
    but WITHOUT ANY WARRANTY; without even the implied warranty of
    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
    Lesser General Public License for more details.

    You should have received a copy of the GNU Lesser General Public License
    (lgpl.txt) along with this library; if not, write to the Free Software
    Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA

    [email protected], [email protected], [email protected],
    [email protected], [email protected], [email protected]

  ***************************************************************** */
/*!
  \file   VertexMover.hpp
  \brief

  The VertexMover Class is the base class for all the smoothing and
  optimizing algorythms

  \author Thomas Leurent
  \date   2002-01-17
*/

#ifndef Mesquite_VertexMover_hpp
#define Mesquite_VertexMover_hpp

#include "Mesquite.hpp"
#include "QualityImprover.hpp"
#include "OFEvaluator.hpp"
#include "MeshInterface.hpp"

namespace MBMesquite
{
class ObjectiveFunction;
class PatchData;
class ParallelMesh;

/*! \class VertexMover
  Base class for all Vertex Movers.
 */
class MESQUITE_EXPORT VertexMover : public QualityImprover
{
  protected:
    VertexMover( ObjectiveFunction* OF = NULL );

  public:
    // virtual destructor ensures use of polymorphism during destruction
    virtual ~VertexMover();

    virtual void initialize_queue( MeshDomainAssoc* mesh_and_domain, const Settings* settings, MsqError& err );<--- Function in derived class

    virtual double loop_over_mesh( MeshDomainAssoc* mesh_and_domain, const Settings* settings, MsqError& err );

    virtual double loop_over_mesh( ParallelMesh* mesh, MeshDomain* domain, const Settings* settings, MsqError& err );

    /**\brief Do Nash-game type optimization
     *
     * Default, opposite of \c do_block_coordinate_descent().
     */
    inline void do_nash_game_optimization()
    {
        objFuncEval.do_nash_game();
    }

    /**\brief Check if optmizer will do Nash-game type optimization
     *
     * Default, opposite of \c is_block_coordinate_descent_optimization().
     */
    inline bool is_nash_game_optimization() const
    {
        return objFuncEval.is_nash_game();
    }

    /**\brief Do block coordinate descent optimization
     *
     * Opposite of \c do_nash_game().
     */
    inline void do_block_coordinate_descent_optimization()
    {
        objFuncEval.do_block_coordinate_descent();
    }

    /**\brief Check if optmizer will do block coordinate descent type optimization
     *
     * Default, opposite of \c is_nash_game_optimization().
     */
    inline bool is_block_coordinate_descent_optimization() const
    {
        return objFuncEval.is_block_coordinate_descent();
    }

    /**\brief Use Jacobi iteration for optimization
     *
     * Opposite of \c do_gauss_optimization()
     */
    inline void do_jacobi_optimization()
    {
        jacobiOpt = true;
    }

    /**\brief Check if optimization will use Jacobi iteration
     *
     * Opposite of \c is_gauss_optimization()
     */
    inline bool is_jacobi_optimization() const
    {
        return jacobiOpt;
    }

    /**\brief Use Gauss-Seidel iteration for optimization
     *
     * Default, opposite of \c do_jacobi_optimization()
     */
    inline void do_gauss_optimization()
    {
        jacobiOpt = false;
    }

    /**\brief Check if optimization will use Gauss-Seidel iteration
     *
     * Default, opposite of \c is_jacobi_optimization()
     */
    inline bool is_gauss_optimization() const
    {
        return !jacobiOpt;
    }

  protected:
    virtual void initialize( PatchData& pd, MsqError& err ) = 0;
    virtual void cleanup()                                  = 0;
    virtual void optimize_vertex_positions( PatchData& pd,
                                            MsqError& err ) = 0;  // modifies the PatchData object

    virtual void initialize_mesh_iteration( PatchData& pd, MsqError& err ) = 0;
    virtual void terminate_mesh_iteration( PatchData&, MsqError& err )     = 0;

    OFEvaluator& get_objective_function_evaluator()
    {
        return objFuncEval;
    }

    static TagHandle get_jacobi_coord_tag( Mesh* mesh, MsqError& err );
    static void commit_jacobi_coords( TagHandle tag, Mesh* mesh, MsqError& err );

  private:
    OFEvaluator objFuncEval;
    bool jacobiOpt;
};

}  // namespace MBMesquite
#endif  // Mesquite_VertexMover_hpp