|
| 1 | +// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*- |
| 2 | +//----------------------------------------------------------------------------- |
| 3 | +// Copyright 2000-2025 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com) |
| 4 | +// See the top-level COPYRIGHT file for details. |
| 5 | +// SPDX-License-Identifier: Apache-2.0 |
| 6 | +//----------------------------------------------------------------------------- |
| 7 | +/*---------------------------------------------------------------------------*/ |
| 8 | +/* PETScDoFLinearSystem.cc (C) 2022-2025 */ |
| 9 | +/* */ |
| 10 | +/* Linear system: Matrix A + Vector x + Vector b for Ax=b. */ |
| 11 | +/*---------------------------------------------------------------------------*/ |
| 12 | +/*---------------------------------------------------------------------------*/ |
| 13 | + |
| 14 | +#include "DoFLinearSystem.h" |
| 15 | + |
| 16 | +#include <mpi.h> |
| 17 | +#include <arccore/trace/TraceAccessor.h> |
| 18 | + |
| 19 | +#include <arcane/accelerator/RunCommandLoop.h> |
| 20 | +#include <arcane/utils/FatalErrorException.h> |
| 21 | +#include <arcane/utils/PlatformUtils.h> |
| 22 | +#include <arcane/utils/ArcaneGlobal.h> |
| 23 | +#include <arcane/utils/MemoryUtils.h> |
| 24 | +#include <arcane/utils/MemoryView.h> |
| 25 | +#include <arcane/utils/ITraceMng.h> |
| 26 | +#include <arcane/utils/NumArray.h> |
| 27 | +#include <arcane/utils/CommandLineArguments.h> |
| 28 | + |
| 29 | +#include <arcane/core/ServiceFactory.h> |
| 30 | +#include <arcane/core/VariableTypes.h> |
| 31 | +#include <arcane/core/BasicService.h> |
| 32 | +#include <arcane/core/IParallelMng.h> |
| 33 | +#include <arcane/core/IItemFamily.h> |
| 34 | +#include <arcane/core/ItemPrinter.h> |
| 35 | +#include <arcane/core/Timer.h> |
| 36 | + |
| 37 | +#include <arcane/accelerator/VariableViews.h> |
| 38 | +#include <arcane/accelerator/core/Runner.h> |
| 39 | +#include <arcane/accelerator/core/Memory.h> |
| 40 | + |
| 41 | +#include "IDoFLinearSystemFactory.h" |
| 42 | +#include "internal/CsrDoFLinearSystemImpl.h" |
| 43 | + |
| 44 | +#include <alien/arcane_tools/accessors/ItemVectorAccessor.h> |
| 45 | +#include <alien/core/block/VBlock.h> |
| 46 | + |
| 47 | +#include <alien/arcane_tools/IIndexManager.h> |
| 48 | +#include <alien/arcane_tools/indexManager/BasicIndexManager.h> |
| 49 | +#include <alien/arcane_tools/indexManager/SimpleAbstractFamily.h> |
| 50 | +#include <alien/arcane_tools/distribution/DistributionFabric.h> |
| 51 | +#include <alien/arcane_tools/indexSet/IndexSetFabric.h> |
| 52 | +#include <alien/arcane_tools/data/Space.h> |
| 53 | + |
| 54 | +#include <alien/kernels/simple_csr/algebra/SimpleCSRLinearAlgebra.h> |
| 55 | +#include <alien/kernels/simple_csr/algebra/SimpleCSRInternalLinearAlgebra.h> |
| 56 | + |
| 57 | +#include <alien/ref/AlienRefSemantic.h> |
| 58 | +#include <alien/ref/AlienImportExport.h> |
| 59 | + |
| 60 | +#include <alien/kernels/redistributor/Redistributor.h> |
| 61 | +#include <alien/ref/data/scalar/RedistributedVector.h> |
| 62 | +#include <alien/ref/data/scalar/RedistributedMatrix.h> |
| 63 | +#include <alien/ref/import_export/MatrixMarketSystemWriter.h> |
| 64 | + |
| 65 | +#include <alien/expression/solver/SolverStater.h> |
| 66 | +#include <alien/AlienLegacyConfig.h> |
| 67 | + |
| 68 | +#include <AlienDoFLinearSystemFactory_axl.h> |
| 69 | + |
| 70 | +#ifdef ALIEN_USE_PETSC |
| 71 | +#include <alien/kernels/petsc/io/AsciiDumper.h> |
| 72 | +#include <alien/kernels/petsc/algebra/PETScLinearAlgebra.h> |
| 73 | +#endif |
| 74 | + |
| 75 | +#ifdef ALIEN_USE_HYPRE |
| 76 | +#include <alien/kernels/hypre/HypreBackEnd.h> |
| 77 | +#include <alien/kernels/hypre/data_structure/HypreMatrix.h> |
| 78 | +#include <alien/kernels/hypre/data_structure/HypreVector.h> |
| 79 | +#include <alien/kernels/hypre/algebra/HypreLinearAlgebra.h> |
| 80 | +#endif |
| 81 | + |
| 82 | +namespace Arcane::FemUtils |
| 83 | +{ |
| 84 | + |
| 85 | +using namespace Arcane; |
| 86 | +using namespace Alien; |
| 87 | + |
| 88 | +class AlienDoFLinearSystemImpl |
| 89 | +: public CsrDoFLinearSystemImpl |
| 90 | +{ |
| 91 | + public: |
| 92 | + |
| 93 | + AlienDoFLinearSystemImpl(IItemFamily* dof_family, const String& solver_name) |
| 94 | + : CsrDoFLinearSystemImpl(dof_family, solver_name) |
| 95 | + , m_dof_matrix_numbering(VariableBuildInfo(dof_family, solver_name + "MatrixNumbering")) |
| 96 | + { |
| 97 | + info() << "[Alien-Info] Creating AlienDoFLinearSystemImpl()"; |
| 98 | + } |
| 99 | + |
| 100 | + ~AlienDoFLinearSystemImpl() override |
| 101 | + { |
| 102 | + info() << "[Alien-Info] Calling AlienDoFLinearSystemImpl destructor"; |
| 103 | + } |
| 104 | + |
| 105 | + public: |
| 106 | + |
| 107 | + void solve() override; |
| 108 | + |
| 109 | + void setSolverCommandLineArguments(const CommandLineArguments& args) override |
| 110 | + { |
| 111 | + info() << "[Alien-Info] initialize command lines arguments"; |
| 112 | + auto argv = *args.commandLineArgv(); |
| 113 | + auto o = info() << "[Alien-Info] ./" << argv[0]; |
| 114 | + |
| 115 | + |
| 116 | + for (int i = 1; i < *args.commandLineArgc(); i++) |
| 117 | + o << ' ' << argv[i]; |
| 118 | + } |
| 119 | + |
| 120 | + void setSolver(Alien::ILinearSolver* s) { m_solver_backend = s; } |
| 121 | + |
| 122 | + CaseOptionsAlienDoFLinearSystemFactory* options; |
| 123 | + |
| 124 | + private: |
| 125 | + |
| 126 | + Alien::ILinearSolver* m_solver_backend; |
| 127 | + VariableDoFInt32 m_dof_matrix_numbering; |
| 128 | + |
| 129 | + NumArray<Real, MDDim1> m_rhs_work_values; |
| 130 | + //! Work array to store values of solution vector in parallel |
| 131 | + NumArray<Real, MDDim1> m_result_work_values; |
| 132 | +}; |
| 133 | + |
| 134 | +void AlienDoFLinearSystemImpl:: |
| 135 | +solve() |
| 136 | +{ |
| 137 | + info() << "[Alien-Info] Calling Alien solver"; |
| 138 | + |
| 139 | + IItemFamily* dof_family = dofFamily(); |
| 140 | + IParallelMng* pm = dof_family->parallelMng(); |
| 141 | + Runner runner = this->runner(); |
| 142 | + CSRFormatView csr_view = this->getCSRValues(); |
| 143 | + |
| 144 | + auto areaU = dof_family->allItems(); |
| 145 | + |
| 146 | + Alien::ArcaneTools::BasicIndexManager index_manager(pm); |
| 147 | + index_manager.setTraceMng(traceMng()); |
| 148 | + |
| 149 | + auto indexSetU = index_manager.buildScalarIndexSet("U", areaU); |
| 150 | + index_manager.prepare(); |
| 151 | + UniqueArray<Arccore::Integer> allUIndex = index_manager.getIndexes(indexSetU); |
| 152 | + |
| 153 | + Alien::ArcaneTools::Space space(&index_manager, "TestSpace"); |
| 154 | + auto mdist = Alien::ArcaneTools::createMatrixDistribution(space); |
| 155 | + auto vdist = Alien::ArcaneTools::createVectorDistribution(space); |
| 156 | + |
| 157 | + Alien::Vector vectorB(vdist); |
| 158 | + Alien::Vector vectorX(vdist); |
| 159 | + Alien::Matrix matrixA(mdist); |
| 160 | + // local matrix for exact measure without side effect |
| 161 | + // (however, you can reuse a matrix with several |
| 162 | + // builder) |
| 163 | + |
| 164 | + Real a1 = platform::getRealTime(); |
| 165 | + |
| 166 | + pm->barrier(); |
| 167 | + |
| 168 | + { |
| 169 | + Alien::MatrixProfiler profiler(matrixA); |
| 170 | + /////////////////////////////////////////////////////////////////////////// |
| 171 | + // |
| 172 | + // DEFINE PROFILE |
| 173 | + // |
| 174 | + ENUMERATE_DOF(idof, dof_family->allItems()) { |
| 175 | + if (!idof->isOwn()) { |
| 176 | + continue; |
| 177 | + } |
| 178 | + |
| 179 | + int i = idof.index(); |
| 180 | + // info() << "owned index: " << i << " local id: " << idof.localId(); |
| 181 | + for (CsrRowColumnIndex csr_index : csr_view.rowRange(i)) { |
| 182 | + Int32 column_index = csr_view.column(csr_index); |
| 183 | + // info() << "column index: " << column_index << " allUIndex: " << allUIndex[idof.localId()] << " csr index: " << csr_index; |
| 184 | + profiler.addMatrixEntry(allUIndex[idof.localId()], allUIndex[column_index]); |
| 185 | + } |
| 186 | + } |
| 187 | + } |
| 188 | + { |
| 189 | + Alien::ProfiledMatrixBuilder builder( |
| 190 | + matrixA, Alien::ProfiledMatrixOptions::eResetValues); |
| 191 | + |
| 192 | + ENUMERATE_DOF(idof, dof_family->allItems()) { |
| 193 | + if (!idof->isOwn()) { |
| 194 | + continue; |
| 195 | + } |
| 196 | + |
| 197 | + int i = idof.index(); |
| 198 | + for (CsrRowColumnIndex csr_index : csr_view.rowRange(i)) { |
| 199 | + Int32 column_index = csr_view.column(csr_index); |
| 200 | + // info() << "row: " << i << " col: " << column_index << " val: " << csr_view.value(csr_index); |
| 201 | + builder(allUIndex[idof.localId()], allUIndex[column_index]) += csr_view.value(csr_index); |
| 202 | + } |
| 203 | + } |
| 204 | + builder.finalize(); |
| 205 | + } |
| 206 | + |
| 207 | + Real a2 = platform::getRealTime(); |
| 208 | + info() << "[Alien-Timer] Time to create matrix = " << (a2 - a1); |
| 209 | + |
| 210 | + VariableDoFReal& rhs_variable = this->rhsVariable(); |
| 211 | + VariableDoFReal& dof_variable = this->solutionVariable(); |
| 212 | + auto rhs_data = rhs_variable.asArray(); |
| 213 | + auto result_data = dof_variable.asArray(); |
| 214 | + |
| 215 | + { |
| 216 | + Alien::VectorWriter writer(vectorX); |
| 217 | + |
| 218 | + ENUMERATE_DOF (idof, areaU.own()) { |
| 219 | + const Integer iIndex = allUIndex[idof->localId()]; |
| 220 | + // info() << "local: " << idof->localId() << " global: " << idof.index(); |
| 221 | + // info() << "local " << idof->localId()<< " iIndex " << iIndex << " res val: " << result_data[idof->localId()]; |
| 222 | + writer[iIndex] = result_data[idof->localId()]; |
| 223 | + } |
| 224 | + } |
| 225 | + |
| 226 | + { |
| 227 | + Alien::VectorWriter writer(vectorB); |
| 228 | + |
| 229 | + ENUMERATE_DOF (idof, areaU.own()) { |
| 230 | + const Integer iIndex = allUIndex[idof->localId()]; |
| 231 | + // info() << "local " << idof->localId( )<< " iIndex " << iIndex << " rhs val: " << rhs_data[idof->localId()]; |
| 232 | + writer[iIndex] = rhs_data[idof->localId()]; |
| 233 | + } |
| 234 | + } |
| 235 | + |
| 236 | + Real a3 = platform::getRealTime(); |
| 237 | + info() << "[Alien-Timer] Time to create vectors = " << (a3 - a2); |
| 238 | + |
| 239 | + m_solver_backend->solve(matrixA, vectorB, vectorX); |
| 240 | + |
| 241 | + Real a4 = platform::getRealTime(); |
| 242 | + info() << "[Alien-Timer] Time to solve = " << (a4 - a3); |
| 243 | + |
| 244 | + Alien::SolverStatus status = m_solver_backend->getStatus(); |
| 245 | + |
| 246 | + if (status.succeeded) { |
| 247 | + info() << "[Alien-Info] " << "Converged in " << status.iteration_count + 1 << " iterations"; |
| 248 | + |
| 249 | + Alien::VectorReader reader(vectorX); |
| 250 | + // TODO understand this |
| 251 | + ENUMERATE_DOF(idof, areaU.own()) { |
| 252 | + const Integer iIndex = allUIndex[idof->localId()]; |
| 253 | + // info() << "val:" << reader[iIndex]; |
| 254 | + dof_variable[idof] = reader[iIndex]; |
| 255 | + } |
| 256 | + } |
| 257 | + else |
| 258 | + info()<<"SOLVER FAILED"; |
| 259 | + m_solver_backend->getSolverStat().print(Universe().traceMng(), status, "Linear Solver : "); |
| 260 | +} |
| 261 | + |
| 262 | +class AlienDoFLinearSystemFactoryService |
| 263 | +: public ArcaneAlienDoFLinearSystemFactoryObject |
| 264 | +{ |
| 265 | + public: |
| 266 | + |
| 267 | + explicit AlienDoFLinearSystemFactoryService(const ServiceBuildInfo& sbi) |
| 268 | + : ArcaneAlienDoFLinearSystemFactoryObject(sbi) |
| 269 | + { |
| 270 | + info() << "[Alien-Info] Create AlienDoF"; |
| 271 | + }; |
| 272 | + IDoFLinearSystemImpl* |
| 273 | + createInstance(ISubDomain* sd, IItemFamily* dof_family, const String& solver_name) override |
| 274 | + { |
| 275 | + auto* x = new AlienDoFLinearSystemImpl(dof_family, solver_name); |
| 276 | + x->options = options(); |
| 277 | + |
| 278 | + Alien::ILinearSolver* solver_backend = options()->linearSolver.instance(); |
| 279 | + x->setSolver(solver_backend); |
| 280 | + solver_backend->init(); |
| 281 | + |
| 282 | + return x; |
| 283 | + } |
| 284 | +}; |
| 285 | + |
| 286 | +ARCANE_REGISTER_SERVICE_ALIENDOFLINEARSYSTEMFACTORY(AlienLinearSystem, |
| 287 | + AlienDoFLinearSystemFactoryService); |
| 288 | + |
| 289 | +} // namespace Arcane::FemUtils |
0 commit comments