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
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
/*
 * compareMaps.cpp
 * this tool will take 2 existing map files in nc format, and compare their sparse matrices
 * we can compare xc, yc, areas, fractions with ncdiff from nco
 * maybe there is another utility in nco, need to ask Charlie Zender
 *
 * example of usage:
 * ./mbcmpmaps -i map1.nc -j map2.nc
 * will look for row, col, S entries, and use eigen3 sparse matrix constructor
 *
 * can be built only if netcdf and eigen3 are available
 *
 *
 */
#include "moab/MOABConfig.h"

#ifndef MOAB_HAVE_EIGEN3
#error compareMaps tool requires eigen3 configuration
#endif

#include "moab/ProgOptions.hpp"

#include "netcdf.h"
#include <cmath>
#include <iomanip>
#include <Eigen/Sparse>
#define ERR_NC( e )                                \
    {                                              \
        printf( "Error: %s\n", nc_strerror( e ) ); \
        exit( 2 );                                 \
    }

// copy from ReadNCDF.cpp some useful macros for reading from a netcdf file (exodus?)
// ncFile is an integer initialized when opening the nc file in read mode

int ncFile1;

#define INS_ID( stringvar, prefix, id ) sprintf( stringvar, prefix, id )

#define GET_DIM1( ncdim, name, val )                            \
    {                                                           \
        int gdfail = nc_inq_dimid( ncFile1, name, &( ncdim ) ); \
        if( NC_NOERR == gdfail )                                \
        {                                                       \
            size_t tmp_val;                                     \
            gdfail = nc_inq_dimlen( ncFile1, ncdim, &tmp_val ); \
            if( NC_NOERR != gdfail )                            \
            {                                                   \
                ERR_NC( gdfail )                                \
            }                                                   \
            else                                                \
                ( val ) = tmp_val;                              \
        }                                                       \
        else                                                    \
            ( val ) = 0;                                        \
    }

#define GET_DIMB1( ncdim, name, varname, id, val ) \
    INS_ID( name, varname, id );                   \
    GET_DIM1( ncdim, name, val );

#define GET_VAR1( name, id, dims )                                     \
    {                                                                  \
        ( id )     = -1;                                               \
        int gvfail = nc_inq_varid( ncFile1, name, &( id ) );           \
        if( NC_NOERR == gvfail )                                       \
        {                                                              \
            int ndims;                                                 \
            gvfail = nc_inq_varndims( ncFile1, id, &ndims );           \
            if( NC_NOERR == gvfail )                                   \
            {                                                          \
                ( dims ).resize( ndims );                              \
                gvfail = nc_inq_vardimid( ncFile1, id, &( dims )[0] ); \
                if( NC_NOERR != gvfail )                               \
                {                                                      \
                    ERR_NC( gvfail )                                   \
                }                                                      \
            }                                                          \
        }                                                              \
    }

#define GET_1D_INT_VAR1( name, id, vals )                                               \
    {                                                                                   \
        GET_VAR1( name, id, vals );                                                     \
        if( -1 != ( id ) )                                                              \
        {                                                                               \
            size_t ntmp;                                                                \
            int ivfail = nc_inq_dimlen( ncFile1, ( vals )[0], &ntmp );                  \
            if( NC_NOERR != ivfail )                                                    \
            {                                                                           \
                ERR_NC( ivfail )                                                        \
            }                                                                           \
            ( vals ).resize( ntmp );                                                    \
            size_t ntmp1 = 0;                                                           \
            ivfail       = nc_get_vara_int( ncFile1, id, &ntmp1, &ntmp, &( vals )[0] ); \
            if( NC_NOERR != ivfail )                                                    \
            {                                                                           \
                ERR_NC( ivfail )                                                        \
            }                                                                           \
        }                                                                               \
    }

#define GET_1D_DBL_VAR1( name, id, vals )                                                  \
    {                                                                                      \
        std::vector< int > dum_dims;                                                       \
        GET_VAR1( name, id, dum_dims );                                                    \
        if( -1 != ( id ) )                                                                 \
        {                                                                                  \
            size_t ntmp;                                                                   \
            int dvfail = nc_inq_dimlen( ncFile1, dum_dims[0], &ntmp );                     \
            if( NC_NOERR != dvfail )                                                       \
            {                                                                              \
                ERR_NC( dvfail )                                                           \
            }                                                                              \
            ( vals ).resize( ntmp );                                                       \
            size_t ntmp1 = 0;                                                              \
            dvfail       = nc_get_vara_double( ncFile1, id, &ntmp1, &ntmp, &( vals )[0] ); \
            if( NC_NOERR != dvfail )                                                       \
            {                                                                              \
                ERR_NC( dvfail )                                                           \
            }                                                                              \
        }                                                                                  \
    }

int ncFile2;

#define GET_DIM2( ncdim, name, val )                            \
    {                                                           \
        int gdfail = nc_inq_dimid( ncFile2, name, &( ncdim ) ); \
        if( NC_NOERR == gdfail )                                \
        {                                                       \
            size_t tmp_val;                                     \
            gdfail = nc_inq_dimlen( ncFile2, ncdim, &tmp_val ); \
            if( NC_NOERR != gdfail )                            \
            {                                                   \
                ERR_NC( gdfail )                                \
            }                                                   \
            else                                                \
                ( val ) = tmp_val;                              \
        }                                                       \
        else                                                    \
            ( val ) = 0;                                        \
    }

#define GET_DIMB2( ncdim, name, varname, id, val ) \
    INS_ID( name, varname, id );                   \
    GET_DIM2( ncdim, name, val );

#define GET_VAR2( name, id, dims )                                     \
    {                                                                  \
        ( id )     = -1;                                               \
        int gvfail = nc_inq_varid( ncFile2, name, &( id ) );           \
        if( NC_NOERR == gvfail )                                       \
        {                                                              \
            int ndims;                                                 \
            gvfail = nc_inq_varndims( ncFile2, id, &ndims );           \
            if( NC_NOERR == gvfail )                                   \
            {                                                          \
                ( dims ).resize( ndims );                              \
                gvfail = nc_inq_vardimid( ncFile2, id, &( dims )[0] ); \
                if( NC_NOERR != gvfail )                               \
                {                                                      \
                    ERR_NC( gvfail )                                   \
                }                                                      \
            }                                                          \
        }                                                              \
    }

#define GET_1D_INT_VAR2( name, id, vals )                                               \
    {                                                                                   \
        GET_VAR2( name, id, vals );                                                     \
        if( -1 != ( id ) )                                                              \
        {                                                                               \
            size_t ntmp;                                                                \
            int ivfail = nc_inq_dimlen( ncFile2, ( vals )[0], &ntmp );                  \
            if( NC_NOERR != ivfail )                                                    \
            {                                                                           \
                ERR_NC( ivfail )                                                        \
            }                                                                           \
            ( vals ).resize( ntmp );                                                    \
            size_t ntmp1 = 0;                                                           \
            ivfail       = nc_get_vara_int( ncFile2, id, &ntmp1, &ntmp, &( vals )[0] ); \
            if( NC_NOERR != ivfail )                                                    \
            {                                                                           \
                ERR_NC( ivfail )                                                        \
            }                                                                           \
        }                                                                               \
    }

#define GET_1D_DBL_VAR2( name, id, vals )                                                  \
    {                                                                                      \
        std::vector< int > dum_dims;                                                       \
        GET_VAR2( name, id, dum_dims );                                                    \
        if( -1 != ( id ) )                                                                 \
        {                                                                                  \
            size_t ntmp;                                                                   \
            int dvfail = nc_inq_dimlen( ncFile2, dum_dims[0], &ntmp );                     \
            if( NC_NOERR != dvfail )                                                       \
            {                                                                              \
                ERR_NC( dvfail )                                                           \
            }                                                                              \
            ( vals ).resize( ntmp );                                                       \
            size_t ntmp1 = 0;                                                              \
            dvfail       = nc_get_vara_double( ncFile2, id, &ntmp1, &ntmp, &( vals )[0] ); \
            if( NC_NOERR != dvfail )                                                       \
            {                                                                              \
                ERR_NC( dvfail )                                                           \
            }                                                                              \
        }                                                                                  \
    }

#define GET_2D_DBL_VAR1( name, id, vals )                                                   \
    {                                                                                       \
        std::vector< int > dum_dims;                                                        \
        GET_VAR1( name, id, dum_dims );                                                     \
        if( -1 != ( id ) )                                                                  \
        {                                                                                   \
            size_t ntmp[2];                                                                 \
            int dvfail = nc_inq_dimlen( ncFile1, dum_dims[0], &ntmp[0] );                   \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
            dvfail = nc_inq_dimlen( ncFile1, dum_dims[1], &ntmp[1] );                       \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
            ( vals ).resize( ntmp[0] * ntmp[1] );                                           \
            size_t ntmp1[2] = { 0, 0 };                                                     \
            dvfail          = nc_get_vara_double( ncFile1, id, ntmp1, ntmp, &( vals )[0] ); \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
        }                                                                                   \
    }

#define GET_2D_DBL_VAR2( name, id, vals )                                                   \
    {                                                                                       \
        std::vector< int > dum_dims;                                                        \
        GET_VAR2( name, id, dum_dims );                                                     \
        if( -1 != ( id ) )                                                                  \
        {                                                                                   \
            size_t ntmp[2];                                                                 \
            int dvfail = nc_inq_dimlen( ncFile2, dum_dims[0], &ntmp[0] );                   \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
            dvfail = nc_inq_dimlen( ncFile2, dum_dims[1], &ntmp[1] );                       \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
            ( vals ).resize( ntmp[0] * ntmp[1] );                                           \
            size_t ntmp1[2] = { 0, 0 };                                                     \
            dvfail          = nc_get_vara_double( ncFile2, id, ntmp1, ntmp, &( vals )[0] ); \
            if( NC_NOERR != dvfail )                                                        \
            {                                                                               \
                ERR_NC( dvfail )                                                            \
            }                                                                               \
        }                                                                                   \
    }

typedef Eigen::Map< Eigen::VectorXd > EigenV;

void diff_vect( const char* var_name, int n )
{
    // compare frac_a between maps
    std::vector< double > fraca1( n ), fraca2( n );
    int idfa1, idfa2;
    GET_1D_DBL_VAR1( var_name, idfa1, fraca1 );
    EigenV fa1( fraca1.data(), n );
    GET_1D_DBL_VAR2( var_name, idfa2, fraca2 );
    EigenV fa2( fraca2.data(), n );

    EigenV diff( fraca2.data(), n );
    diff = fa1 - fa2;

    int imin, imax;
    double minV = diff.minCoeff( &imin );
    double maxV = diff.maxCoeff( &imax );
    std::cout << var_name << " diff norm: " << diff.norm() << " min at " << imin << " : " << minV << " max at " << imax
              << " : " << maxV << "\n";
    return;
}

void diff_2d_vect( const char* var_name, int n )
{
    // compare frac_a between maps
    std::vector< double > fraca1( n ), fraca2( n );
    int idfa1, idfa2;
    GET_2D_DBL_VAR1( var_name, idfa1, fraca1 );
    EigenV fa1( fraca1.data(), n );
    GET_2D_DBL_VAR2( var_name, idfa2, fraca2 );
    EigenV fa2( fraca2.data(), n );
    std::cout << var_name << " diff norm: " << ( fa1 - fa2 ).norm() << "\n";
    return;
}
int main( int argc, char* argv[] )
{

    ProgOptions opts;

    std::string inputfile1, inputfile2;
    opts.addOpt< std::string >( "firstMap,i", "input filename 1", &inputfile1 );
    opts.addOpt< std::string >( "secondMap,j", "input second map", &inputfile2 );
    int print_diff = 20;
    opts.addOpt< int >( "print_differences,p", "print differences ", &print_diff );
    double fraction = 0.9;
    opts.addOpt< double >( "fraction_diff,f", "fraction threshold", &fraction );

    opts.parseCommandLine( argc, argv );

    // Open netcdf/exodus file
    int fail = nc_open( inputfile1.c_str(), 0, &ncFile1 );
    if( NC_NOWRITE != fail )
    {
        ERR_NC( fail )
    }

    std::cout << " opened " << inputfile1 << " for map 1 \n";

    fail = nc_open( inputfile2.c_str(), 0, &ncFile2 );
    if( NC_NOWRITE != fail )
    {
        ERR_NC( fail )
    }

    std::cout << " opened " << inputfile2 << " for map 2 \n";
    int temp_dim;
    int na1, nb1, ns1, na2, nb2, ns2, nv_a, nv_b, nv_a2, nv_b2;
    GET_DIM1( temp_dim, "n_a", na1 );
    GET_DIM2( temp_dim, "n_a", na2 );
    GET_DIM1( temp_dim, "n_b", nb1 );
    GET_DIM2( temp_dim, "n_b", nb2 );
    GET_DIM1( temp_dim, "n_s", ns1 );
    GET_DIM2( temp_dim, "n_s", ns2 );
    GET_DIM1( temp_dim, "nv_a", nv_a );
    GET_DIM2( temp_dim, "nv_a", nv_a2 );
    GET_DIM1( temp_dim, "nv_b", nv_b );
    GET_DIM2( temp_dim, "nv_b", nv_b2 );
    if( nv_a != nv_a2 || nv_b != nv_b2 )
    {
        std::cout << " different nv dimensions:" << nv_a << " == " << nv_a2 << "  or " << nv_b << " == " << nv_b2
                  << "  bail out \n";
        return 1;
    }
    std::cout << " n_a, n_b, n_s : " << na1 << ", " << nb1 << ", " << ns1 << " for map 1 \n";

    std::cout << " n_a, n_b, n_s : " << na2 << ", " << nb2 << ", " << ns2 << " for map 2 \n";
    if( na1 != na2 || nb1 != nb2 )
    {
        std::cout << " different dimensions bail out \n";
        return 1;
    }
    std::vector< int > col1( ns1 ), row1( ns1 );
    std::vector< int > col2( ns2 ), row2( ns2 );

    std::vector< double > val1( ns1 ), val2( ns2 );

    int idrow1, idcol1, idrow2, idcol2, ids1, ids2;
    GET_1D_INT_VAR1( "row", idrow1, row1 );
    GET_1D_INT_VAR1( "col", idcol1, col1 );
    GET_1D_DBL_VAR1( "S", ids1, val1 );
    GET_1D_INT_VAR2( "row", idrow2, row2 );
    GET_1D_INT_VAR2( "col", idcol2, col2 );
    GET_1D_DBL_VAR2( "S", ids2, val2 );

    // first matrix
    typedef Eigen::Triplet< double > Triplet;
    std::vector< Triplet > tripletList;
    tripletList.reserve( ns1 );
    for( int iv = 0; iv < ns1; iv++ )
    {
        tripletList.push_back( Triplet( row1[iv] - 1, col1[iv] - 1, val1[iv] ) );
    }
    Eigen::SparseMatrix< double > weight1( nb1, na1 );

    weight1.setFromTriplets( tripletList.begin(), tripletList.end() );
    weight1.makeCompressed();

    if( ns1 != ns2 ) tripletList.resize( ns2 );
    for( int iv = 0; iv < ns2; iv++ )
    {
        tripletList[iv] = Triplet( row2[iv] - 1, col2[iv] - 1, val2[iv] );
    }
    Eigen::SparseMatrix< double > weight2( nb2, na2 );

    weight2.setFromTriplets( tripletList.begin(), tripletList.end() );
    weight2.makeCompressed();

    // default storage type is column major
    Eigen::SparseMatrix< double > diff = weight1 - weight2;
    diff.makeCompressed();  // is it needed or not ?
    auto coeffs = diff.coeffs();
    double maxv = coeffs.maxCoeff();
    double minv = coeffs.minCoeff();
    std::cout << " euclidian norm for difference: " << diff.norm()
              << " \n squared norm for difference: " << diff.squaredNorm() << "\n"
              << " minv: " << minv << " maxv: " << maxv << "\n";
    // print out the first 10 positions for which the value is outside 90% of min/max values
    double min_threshold = fraction * minv;
    double max_threshold = fraction * maxv;
    int counter          = 0;
    std::cout << std::setprecision( 9 );
    for( int k = 0; ( k < diff.outerSize() ); ++k )  // this is by column
    {
        for( Eigen::SparseMatrix< double >::InnerIterator it( diff, k ); ( it ) && ( counter < print_diff ); ++it )
        {
            double val = it.value();
            if( val <= min_threshold || val >= max_threshold )
            {
                int row = it.row();
                int col = it.col();
                std::cout << " counter:" << counter << "\t col: " << col + 1 << "\t row: " << row + 1
                          << "\t diff: " << val;
                std::cout << "\t map1: " << weight1.coeffRef( row, col ) << "\t map2: " << weight2.coeffRef( row, col )
                          << "\n";  // row index
                counter++;
            }
        }
    }

    // compare frac_a between maps
    diff_vect( "frac_a", na1 );
    diff_vect( "frac_b", nb1 );
    diff_vect( "area_a", na1 );
    diff_vect( "area_b", nb1 );
    diff_vect( "yc_a", na1 );
    diff_vect( "yc_b", nb1 );
    diff_vect( "xc_a", na1 );
    diff_vect( "xc_b", nb1 );
    diff_2d_vect( "xv_a", na1 * nv_a );
    diff_2d_vect( "yv_a", na1 * nv_a );
    diff_2d_vect( "xv_b", nb1 * nv_b );
    diff_2d_vect( "yv_b", nb1 * nv_b );

    return 0;
}