// Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced // at the Lawrence Livermore National Laboratory. All Rights reserved. See files // LICENSE and NOTICE for details. LLNL-CODE-806117. // // This file is part of the MFEM library. For more information and source code // availability visit https://mfem.org. // // MFEM is free software; you can redistribute it and/or modify it under the // terms of the BSD-3 license. We welcome feedback and contributions, see file // CONTRIBUTING.md for details. // Implementations of classes FABilinearFormExtension, EABilinearFormExtension, // PABilinearFormExtension and MFBilinearFormExtension. #include "nonlinearform.hpp" #include "ceed/interface/util.hpp" namespace mfem { NonlinearFormExtension::NonlinearFormExtension(const NonlinearForm *nlf) : Operator(nlf->FESpace()->GetVSize()), nlf(nlf) { } PANonlinearFormExtension::PANonlinearFormExtension(const NonlinearForm *nlf): NonlinearFormExtension(nlf), fes(*nlf->FESpace()), dnfi(*nlf->GetDNFI()), elemR(nullptr), Grad(*this) { if (!DeviceCanUseCeed()) { elemR = fes.GetElementRestriction(ElementDofOrdering::LEXICOGRAPHIC); // TODO: optimize for the case when 'elemR' is identity xe.SetSize(elemR->Height(), Device::GetMemoryType()); ye.SetSize(elemR->Height(), Device::GetMemoryType()); } ye.UseDevice(true); } double PANonlinearFormExtension::GetGridFunctionEnergy(const Vector &x) const { double energy = 0.0; elemR->Mult(x, xe); for (int i = 0; i < dnfi.Size(); i++) { energy += dnfi[i]->GetLocalStateEnergyPA(xe); } return energy; } void PANonlinearFormExtension::Assemble() { MFEM_VERIFY(nlf->GetInteriorFaceIntegrators().Size() == 0 && nlf->GetBdrFaceIntegrators().Size() == 0, "face integrators are not supported yet"); for (int i = 0; i < dnfi.Size(); ++i) { dnfi[i]->AssemblePA(fes); } } void PANonlinearFormExtension::Mult(const Vector &x, Vector &y) const { if (!DeviceCanUseCeed()) { ye = 0.0; elemR->Mult(x, xe); for (int i = 0; i < dnfi.Size(); ++i) { dnfi[i]->AddMultPA(xe, ye); } elemR->MultTranspose(ye, y); } else { y.UseDevice(true); // typically this is a large vector, so store on device y = 0.0; for (int i = 0; i < dnfi.Size(); ++i) { dnfi[i]->AddMultPA(x, y); } } } Operator &PANonlinearFormExtension::GetGradient(const Vector &x) const { Grad.AssembleGrad(x); return Grad; } void PANonlinearFormExtension::Update() { height = width = fes.GetVSize(); elemR = fes.GetElementRestriction(ElementDofOrdering::LEXICOGRAPHIC); xe.SetSize(elemR->Height()); ye.SetSize(elemR->Height()); Grad.Update(); } PANonlinearFormExtension::Gradient::Gradient(const PANonlinearFormExtension &e): Operator(e.Height()), ext(e) { } void PANonlinearFormExtension::Gradient::AssembleGrad(const Vector &g) { ext.elemR->Mult(g, ext.xe); for (int i = 0; i < ext.dnfi.Size(); ++i) { ext.dnfi[i]->AssembleGradPA(ext.xe, ext.fes); } } void PANonlinearFormExtension::Gradient::Mult(const Vector &x, Vector &y) const { ext.ye = 0.0; ext.elemR->Mult(x, ext.xe); for (int i = 0; i < ext.dnfi.Size(); ++i) { ext.dnfi[i]->AddMultGradPA(ext.xe, ext.ye); } ext.elemR->MultTranspose(ext.ye, y); } void PANonlinearFormExtension::Gradient::AssembleDiagonal(Vector &diag) const { MFEM_ASSERT(diag.Size() == Height(), "Vector for holding diagonal has wrong size!"); ext.ye = 0.0; for (int i = 0; i < ext.dnfi.Size(); ++i) { ext.dnfi[i]->AssembleGradDiagonalPA(ext.ye); } ext.elemR->MultTranspose(ext.ye, diag); } void PANonlinearFormExtension::Gradient::Update() { height = width = ext.Height(); } MFNonlinearFormExtension::MFNonlinearFormExtension(const NonlinearForm *form): NonlinearFormExtension(form), fes(*form->FESpace()) { if (!DeviceCanUseCeed()) { const ElementDofOrdering ordering = ElementDofOrdering::LEXICOGRAPHIC; elem_restrict_lex = fes.GetElementRestriction(ordering); if (elem_restrict_lex) // replace with a check for not identity { localX.SetSize(elem_restrict_lex->Height(), Device::GetMemoryType()); localY.SetSize(elem_restrict_lex->Height(), Device::GetMemoryType()); localY.UseDevice(true); // ensure 'localY = 0.0' is done on device } } } void MFNonlinearFormExtension::Assemble() { const Array &integrators = *nlf->GetDNFI(); const int Ni = integrators.Size(); for (int i = 0; i < Ni; ++i) { integrators[i]->AssembleMF(fes); } } void MFNonlinearFormExtension::Mult(const Vector &x, Vector &y) const { const Array &integrators = *nlf->GetDNFI(); const int iSz = integrators.Size(); // replace the check 'elem_restrict_lex' with a check for not identity if (elem_restrict_lex && !DeviceCanUseCeed()) { elem_restrict_lex->Mult(x, localX); localY = 0.0; for (int i = 0; i < iSz; ++i) { integrators[i]->AddMultMF(localX, localY); } elem_restrict_lex->MultTranspose(localY, y); } else { y.UseDevice(true); // typically this is a large vector, so store on device y = 0.0; for (int i = 0; i < iSz; ++i) { integrators[i]->AddMultMF(x, y); } } } void MFNonlinearFormExtension::Update() { height = width = fes.GetVSize(); const ElementDofOrdering ordering = ElementDofOrdering::LEXICOGRAPHIC; elem_restrict_lex = fes.GetElementRestriction(ordering); if (elem_restrict_lex) // replace with a check for not identity { localX.SetSize(elem_restrict_lex->Height(), Device::GetMemoryType()); localY.SetSize(elem_restrict_lex->Height(), Device::GetMemoryType()); } } } // namespace mfem