File: LocalAdaptCellRegression.cpp

package info (click to toggle)
stopt 5.12%2Bdfsg-3
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 8,860 kB
  • sloc: cpp: 70,456; python: 5,950; makefile: 72; sh: 57
file content (78 lines) | stat: -rw-r--r-- 2,711 bytes parent folder | download | duplicates (3)
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
// Copyright (C) 2019 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#include "StOpt/regression/LocalAdaptCellRegression.h"

using namespace std;
using namespace Eigen;

namespace StOpt
{
LocalAdaptCellRegression::LocalAdaptCellRegression(const ArrayXi   &p_nbMesh,
        const bool &p_bRotationAndRecale):
    BaseRegression(false, p_bRotationAndRecale),
    m_nbMesh(p_nbMesh),
    m_nbMeshTotal(m_nbMesh.prod()) {}

LocalAdaptCellRegression::LocalAdaptCellRegression(const bool &p_bZeroDate,
        const ArrayXXd  &p_particles,
        const ArrayXi   &p_nbMesh,
        const bool &p_bRotationAndRecale):
    BaseRegression(p_bZeroDate, p_particles, p_bRotationAndRecale),
    m_nbMesh(p_nbMesh),
    m_nbMeshTotal(m_nbMesh.prod()),
    m_simToCell(p_particles.cols())
{}

LocalAdaptCellRegression::LocalAdaptCellRegression(const bool &p_bZeroDate,
        const  ArrayXi &p_nbMesh,
        const  Array< array< double, 2>, Dynamic, Dynamic >   &p_mesh,
        const   ArrayXd &p_meanX,
        const   ArrayXd   &p_etypX, const   MatrixXd   &p_svdMatrix, const  bool &p_bRotationAndRecale) :
    BaseRegression(p_bZeroDate, p_meanX,  p_etypX,  p_svdMatrix, p_bRotationAndRecale), m_nbMesh(p_nbMesh), m_mesh(p_mesh)
{
    if ((!m_bZeroDate) && (p_nbMesh.size() > 0))
    {
        m_nbMeshTotal = m_nbMesh.prod();

        if (!p_bRotationAndRecale)
        {
            m_meanX = ArrayXd::Zero(p_nbMesh.size());
            m_etypX = ArrayXd::Constant(p_nbMesh.size(), 1.);
            m_svdMatrix = MatrixXd::Identity(p_nbMesh.size(), p_nbMesh.size());
        }
    }
    else
        m_nbMeshTotal = 1;
}

LocalAdaptCellRegression::LocalAdaptCellRegression(const LocalAdaptCellRegression   &p_object) : BaseRegression(p_object), m_nbMesh(p_object.getNbMesh()),
    m_nbMeshTotal(p_object.getNbMeshTotal()), m_mesh(p_object.getMesh()),
    m_simToCell(p_object.getSimToCell()),
    m_simulBelongingToCell(p_object.getSimulBelongingToCell())
{
}


void  LocalAdaptCellRegression::evaluateSimulBelongingToCell()
{
    int nbCells = getNumberOfCells();
    m_simulBelongingToCell.resize(nbCells);
    if (m_particles.cols() > 0)
    {
        for (int icell = 0; icell < nbCells; ++icell)
        {
            m_simulBelongingToCell[icell] = make_shared< vector< int> >();
            m_simulBelongingToCell[icell]->reserve(2 * m_particles.cols() / nbCells);
        }
        for (int is = 0; is <  m_simToCell.size(); ++is)
            m_simulBelongingToCell[m_simToCell(is)]->push_back(is);
    }
    else
    {
        m_simulBelongingToCell[0] = make_shared< vector< int> >(1);
        (*m_simulBelongingToCell[0])[0] = 0;
    }
}

}