Skip to content

Commit

Permalink
Add VectorTimeKernel and VectorTimeDerivative
Browse files Browse the repository at this point in the history
  • Loading branch information
cticenhour committed Aug 26, 2019
1 parent db81c8f commit 07c8a8e
Show file tree
Hide file tree
Showing 4 changed files with 189 additions and 0 deletions.
32 changes: 32 additions & 0 deletions framework/include/kernels/VectorTimeDerivative.h
@@ -0,0 +1,32 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "VectorTimeKernel.h"

// Forward Declaration
class VectorTimeDerivative;

template <>
InputParameters validParams<VectorTimeDerivative>();

class VectorTimeDerivative : public VectorTimeKernel
{
public:
VectorTimeDerivative(const InputParameters & parameters);

virtual void computeJacobian() override;

protected:
virtual Real computeQpResidual() override;
virtual Real computeQpJacobian() override;

bool _lumping;
};
37 changes: 37 additions & 0 deletions framework/include/kernels/VectorTimeKernel.h
@@ -0,0 +1,37 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "VectorKernel.h"

// Forward Declaration
class VectorTimeKernel;

template <>
InputParameters validParams<VectorTimeKernel>();

/**
* All vector time kernels should inherit from this class
*
*/
class VectorTimeKernel : public VectorKernel
{
public:
VectorTimeKernel(const InputParameters & parameters);

virtual void computeResidual() override;

protected:
/// Time derivative of u
const VectorVariableValue & _u_dot;

/// Derivative of u_dot with respect to u
const VariableValue & _du_dot_du;
};
66 changes: 66 additions & 0 deletions framework/src/kernels/VectorTimeDerivative.C
@@ -0,0 +1,66 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "VectorTimeDerivative.h"

// MOOSE includes
#include "Assembly.h"
#include "MooseVariableFE.h"

#include "libmesh/quadrature.h"

registerMooseObject("MooseApp", VectorTimeDerivative);

template <>
InputParameters
validParams<VectorTimeDerivative>()
{
InputParameters params = validParams<VectorTimeKernel>();
params.addClassDescription("The time derivative operator with the weak form of $(\\vec{\\psi_i}, "
"\\frac{\\partial \\vec{u_h}}{\\partial t})$.");
params.addParam<bool>("lumping", false, "True for mass matrix lumping, false otherwise");

return params;
}

VectorTimeDerivative::VectorTimeDerivative(const InputParameters & parameters)
: VectorTimeKernel(parameters), _lumping(getParam<bool>("lumping"))
{
}

Real
VectorTimeDerivative::computeQpResidual()
{
return _test[_i][_qp] * _u_dot[_qp];
}

Real
VectorTimeDerivative::computeQpJacobian()
{
return _test[_i][_qp] * _phi[_j][_qp] * _du_dot_du[_qp];
}

void
VectorTimeDerivative::computeJacobian()
{
if (_lumping)
{
prepareMatrixTag(_assembly, _var.number(), _var.number());

precalculateJacobian();
for (_i = 0; _i < _test.size(); _i++)
for (_j = 0; _j < _phi.size(); _j++)
for (_qp = 0; _qp < _qrule->n_points(); _qp++)
_local_ke(_i, _i) += _JxW[_qp] * _coord[_qp] * computeQpJacobian();

accumulateTaggedLocalMatrix();
}
else
VectorTimeKernel::computeJacobian();
}
54 changes: 54 additions & 0 deletions framework/src/kernels/VectorTimeKernel.C
@@ -0,0 +1,54 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "VectorTimeKernel.h"

// MOOSE includes
#include "Assembly.h"
#include "MooseVariableFE.h"
#include "SystemBase.h"

#include "libmesh/quadrature.h"

template <>
InputParameters
validParams<VectorTimeKernel>()
{
InputParameters params = validParams<VectorKernel>();

params.set<MultiMooseEnum>("vector_tags") = "time";
params.set<MultiMooseEnum>("matrix_tags") = "system time";

return params;
}

VectorTimeKernel::VectorTimeKernel(const InputParameters & parameters)
: VectorKernel(parameters), _u_dot(_var.uDot()), _du_dot_du(_var.duDotDu())
{
}

void
VectorTimeKernel::computeResidual()
{
prepareVectorTag(_assembly, _var.number());

precalculateResidual();
for (_i = 0; _i < _test.size(); _i++)
for (_qp = 0; _qp < _qrule->n_points(); _qp++)
_local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual();

accumulateTaggedLocalResidual();

if (_has_save_in)
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
for (unsigned int i = 0; i < _save_in.size(); i++)
_save_in[i]->sys().solution().add_vector(_local_re, _save_in[i]->dofIndices());
}
}

0 comments on commit 07c8a8e

Please sign in to comment.