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
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
/* *****************************************************************
    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]

  ***************************************************************** */
// -*- Mode : c++; tab-width: 3; c-tab-always-indent: t; indent-tabs-mode: nil; c-basic-offset: 3
// -*-
//
//   SUMMARY:
//     USAGE:
//
//    AUTHOR: Michael Brewer
//       ORG: Sandia National Labs
//    E-MAIL: [email protected]
//
// ORIG-DATE: Jan. 29, 2003
//  LAST-MOD: 25-Feb-04 at 10:49:04 by Thomas Leurent
//
// DESCRIPTION:
// ============
/*! \file SphericalGeometryTest.cpp

Regression testing using the spherical geometry capabilities in
SimplifiedGeometryEngine.
 */
// DESCRIP-END.
//

#include "meshfiles.h"

#include "PatchDataInstances.hpp"
#include "cppunit/extensions/HelperMacros.h"
#include <cmath>

#include "Mesquite.hpp"
#include "MsqError.hpp"
#include "Vector3D.hpp"
#include "InstructionQueue.hpp"
#include "PatchData.hpp"
//#include "StoppingCriterion.hpp"
#include "QualityAssessor.hpp"

#include "IdealWeightInverseMeanRatio.hpp"
#include "ConditionNumberQualityMetric.hpp"
#include "LPtoPTemplate.hpp"
#include "EdgeLengthQualityMetric.hpp"
#include "LaplacianSmoother.hpp"
#include "SmartLaplacianSmoother.hpp"
#include "LInfTemplate.hpp"
#include "SteepestDescent.hpp"
#include "ConjugateGradient.hpp"
#include "AspectRatioGammaQualityMetric.hpp"
#include "UntangleBetaQualityMetric.hpp"
#include "MultiplyQualityMetric.hpp"
#include "SphericalDomain.hpp"

#include "MeshImpl.hpp"

#include <iostream>
using std::cout;
using std::endl;
using namespace MBMesquite;

class SphericalGeometryTest : public CppUnit::TestFixture
{
  private:
    CPPUNIT_TEST_SUITE( SphericalGeometryTest );
    // run cg on the quad sphere mesh with condition number L2
    CPPUNIT_TEST( test_cg_mesh_cond_sphere );
    // run cg on the quad sphere mesh with condition number L2
    CPPUNIT_TEST( test_smart_lapl_sphere );
    // run laplacian smoothing on the geo tri mesh
    CPPUNIT_TEST( test_lapl_geo_sphere );

    CPPUNIT_TEST_SUITE_END();

  private:
    double qualTol;  // double used for double comparisons
    int pF;          // PRINT_FLAG
  public:
    void setUp()
    {
        // pF=1;//PRINT_FLAG IS ON
        pF = 0;  // PRINT_FLAG IS OFF
                 // tolerance double
        qualTol = MSQ_MIN;
    }

    void tearDown() {}

  public:
    SphericalGeometryTest() {}<--- Member variable 'SphericalGeometryTest::qualTol' is not initialized in the constructor.<--- Member variable 'SphericalGeometryTest::pF' is not initialized in the constructor.

    void test_cg_mesh_cond_sphere()
    {
        MBMesquite::MeshImpl mesh;
        MBMesquite::MsqPrintError err( cout );
        mesh.read_vtk( MESH_FILES_DIR "2D/vtk/quads/untangled/quads_on_sphere_529.vtk", err );
        CPPUNIT_ASSERT( !err );

        // create geometry: sphere, center (2,2,0), radius 3
        Vector3D center( 2, 2, 0 );
        SphericalDomain msq_geom( center, 3.0 );

        // creates an intruction queue
        InstructionQueue queue1;

        // creates a mean ratio quality metric ...
        ConditionNumberQualityMetric shape;
        UntangleBetaQualityMetric untan;

        // ... and builds an objective function with it
        LPtoPTemplate obj_func( &shape, 2, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // creates the steepest descent optimization procedures
        ConjugateGradient pass1( &obj_func, err );
        // SteepestDescent* pass2 = new SteepestDescent( obj_func );
        pass1.use_global_patch();
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        QualityAssessor qa = QualityAssessor( &shape );

        //**********Set stopping criterion  5 iterates ****************
        // StoppingCriterion sc5(StoppingCriterion::NUMBER_OF_PASSES,5);
        // pass1->set_stopping_criterion(&sc5);
        TerminationCriterion sc5;
        sc5.add_iteration_limit( 5 );
        pass1.set_inner_termination_criterion( &sc5 );
        // CG's debugging print, increase integer to get more print info
        pass1.set_debugging_level( 0 );

        //  queue1.add_preconditioner(pass2, err); CPPUNIT_ASSERT(!err);
        queue1.set_master_quality_improver( &pass1, err );
        CPPUNIT_ASSERT( !err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // launches optimization on mesh_set1
        MeshDomainAssoc mesh_and_domain = MeshDomainAssoc( &mesh, &msq_geom );
        double orig_qa_val              = qa.loop_over_mesh( &mesh_and_domain, 0, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        queue1.run_instructions( &mesh_and_domain, err );
        CPPUNIT_ASSERT( !err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        double fin_qa_val = qa.loop_over_mesh( &mesh_and_domain, 0, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // make sure 'quality' improved
        CPPUNIT_ASSERT( fin_qa_val <= orig_qa_val );
    }
    void test_smart_lapl_sphere()
    {
        MBMesquite::MeshImpl mesh;
        MBMesquite::MsqPrintError err( cout );
        mesh.read_vtk( MESH_FILES_DIR "2D/vtk/quads/untangled/quads_on_sphere_529.vtk", err );

        // create geometry sphere:  ratius 1, centered at (0,0,0)
        Vector3D center( 2, 2, 0 );
        SphericalDomain msq_geom( center, 3.0 );

        // creates an intruction queue
        InstructionQueue queue1;

        // creates an edge length metric ...
        IdealWeightInverseMeanRatio shape_metric( err );
        LInfTemplate shape_func( &shape_metric );

        // create the smart laplacian smoother
        SmartLaplacianSmoother s_lapl( &shape_func );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );

        //*******Set stopping criterion 5 iterates  ***********
        TerminationCriterion sc5;
        sc5.add_iteration_limit( 5 );
        s_lapl.set_outer_termination_criterion( &sc5 );
        // qa, qi, qa
        queue1.set_master_quality_improver( &s_lapl, err );
        CPPUNIT_ASSERT( !err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // launches optimization on mesh_set1
        QualityAssessor qa              = QualityAssessor( &shape_metric );
        MeshDomainAssoc mesh_and_domain = MeshDomainAssoc( &mesh, &msq_geom );
        double orig_val                 = qa.loop_over_mesh( &mesh_and_domain, 0, err );

        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        queue1.run_instructions( &mesh_and_domain, err );
        CPPUNIT_ASSERT( !err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );

        double final_val = qa.loop_over_mesh( &mesh_and_domain, 0, err );

        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // make sure 'quality' improved
        CPPUNIT_ASSERT( final_val < orig_val );
    }

    void test_lapl_geo_sphere()
    {
        MBMesquite::MeshImpl mesh;
        MBMesquite::MsqPrintError err( cout );

        mesh.read_vtk( MESH_FILES_DIR "2D/vtk/tris/untangled/Mesquite_geo_10242.vtk", err );

        // create geometry sphere:  ratius 1, centered at (0,0,0)
        Vector3D center( 0, 0, 0 );
        MBMesquite::SphericalDomain msq_geom( center, 1.0 );

        // creates an intruction queue
        InstructionQueue queue1;

        // creates an edge length metric ...
        EdgeLengthQualityMetric edg_len;

        // create the laplacian smoother
        LaplacianSmoother lapl;
        // Make sure no errors
        CPPUNIT_ASSERT( !err );

        // create a quality assessor
        QualityAssessor qa = QualityAssessor( &edg_len );

        //*******Set stopping criterion 10 iterates  ***********
        // StoppingCriterion sc10(StoppingCriterion::NUMBER_OF_PASSES,10);
        // lapl->set_stopping_criterion(&sc10);
        TerminationCriterion sc10;
        sc10.add_iteration_limit( 10 );
        lapl.set_outer_termination_criterion( &sc10 );
        // qa, qi, qa
        queue1.set_master_quality_improver( &lapl, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // launches optimization on mesh_set1
        MeshDomainAssoc mesh_and_domain = MeshDomainAssoc( &mesh, &msq_geom );
        double orig_qa_val              = qa.loop_over_mesh( &mesh_and_domain, 0, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        queue1.run_instructions( &mesh_and_domain, err );
        CPPUNIT_ASSERT( !err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        double fin_qa_val = qa.loop_over_mesh( &mesh_and_domain, 0, err );
        // Make sure no errors
        CPPUNIT_ASSERT( !err );
        // make sure 'quality' improved
        CPPUNIT_ASSERT( fin_qa_val <= orig_qa_val );
    }
};

CPPUNIT_TEST_SUITE_NAMED_REGISTRATION( SphericalGeometryTest, "SphericalGeometryTest" );
CPPUNIT_TEST_SUITE_NAMED_REGISTRATION( SphericalGeometryTest, "Regression" );