-
Notifications
You must be signed in to change notification settings - Fork 1k
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Add test objects for solving the biharmonic problem.
These Kernels/BCs are needed for properly testing the LaplacianJumpIndicator, but they are probably general enough that they could be added to MOOSE proper, if desired. Refs #2190.
- Loading branch information
1 parent
f599233
commit 9d98df6
Showing
7 changed files
with
277 additions
and
0 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,44 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
|
||
#ifndef BIHARMONICLAPBC_H | ||
#define BIHARMONICLAPBC_H | ||
|
||
#include "IntegratedBC.h" | ||
|
||
// Forward Declarations | ||
class BiharmonicLapBC; | ||
|
||
template <> | ||
InputParameters validParams<BiharmonicLapBC>(); | ||
|
||
/** | ||
* The weak form of the biharmonic equation has a term | ||
* \int -Lap(u) * dv/dn ds | ||
* which we use to weakly impose the value of Lap(u) on the boundary. | ||
*/ | ||
class BiharmonicLapBC : public IntegratedBC | ||
{ | ||
public: | ||
BiharmonicLapBC(const InputParameters & parameters); | ||
virtual ~BiharmonicLapBC() {} | ||
|
||
protected: | ||
virtual Real computeQpResidual(); | ||
|
||
/// User-provided function which computes the Laplacian. | ||
Function & _lap_u; | ||
}; | ||
|
||
#endif |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,57 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
#ifndef FUNCTIONPENALTYFLUXBC_H | ||
#define FUNCTIONPENALTYFLUXBC_H | ||
|
||
#include "IntegratedBC.h" | ||
|
||
class FunctionPenaltyFluxBC; | ||
class Function; | ||
|
||
template <> | ||
InputParameters validParams<FunctionPenaltyFluxBC>(); | ||
|
||
/** | ||
* Penalizes the difference between the current flux and desired flux, | ||
* similarly to penalty Dirichlet BCs. Implements the term: | ||
* | ||
* \int (p * (du/dn - dg/dn) * dv/dn) dx | ||
* | ||
* where p is the (large) penalty parameter, du/dn is the normal flux, | ||
* and dg/dn is the normal component of the gradient of the true solution. | ||
* | ||
* We allow the user to provide the components of the true flux, and | ||
* then compute g for them by dotting those components with the | ||
* outward unit normal. This class is designed to impose a given value | ||
* of the flux as an essential BC in the biharmonic problem. | ||
*/ | ||
class FunctionPenaltyFluxBC : public IntegratedBC | ||
{ | ||
public: | ||
/** | ||
* Factory constructor, takes parameters so that all derived classes can be built using the same | ||
* constructor. | ||
*/ | ||
FunctionPenaltyFluxBC(const InputParameters & parameters); | ||
|
||
protected: | ||
virtual Real computeQpResidual() override; | ||
virtual Real computeQpJacobian() override; | ||
|
||
private: | ||
Function & _func; | ||
Real _p; | ||
}; | ||
|
||
#endif |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,45 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
#ifndef BIHARMONIC_H | ||
#define BIHARMONIC_H | ||
|
||
#include "Kernel.h" | ||
|
||
// Forward Declarations | ||
class Biharmonic; | ||
|
||
template <> | ||
InputParameters validParams<Biharmonic>(); | ||
|
||
/** | ||
* Computes the residual and Jacobian contribution for the weak form | ||
* of the biharmonic equation: | ||
* | ||
* \int Laplacian(u) * Laplacian(v) dx | ||
*/ | ||
class Biharmonic : public Kernel | ||
{ | ||
public: | ||
Biharmonic(const InputParameters & parameters); | ||
|
||
protected: | ||
virtual Real computeQpResidual(); | ||
virtual Real computeQpJacobian(); | ||
|
||
const VariableSecond & _second_u; | ||
const VariablePhiSecond & _second_phi; | ||
const VariableTestSecond & _second_test; | ||
}; | ||
|
||
#endif |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,37 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
|
||
#include "BiharmonicLapBC.h" | ||
#include "Function.h" | ||
|
||
template <> | ||
InputParameters | ||
validParams<BiharmonicLapBC>() | ||
{ | ||
InputParameters params = validParams<IntegratedBC>(); | ||
params.addParam<FunctionName>( | ||
"laplacian_function", "0", "A function representing the weakly-imposed Laplacian."); | ||
return params; | ||
} | ||
|
||
BiharmonicLapBC::BiharmonicLapBC(const InputParameters & parameters) | ||
: IntegratedBC(parameters), _lap_u(getFunction("laplacian_function")) | ||
{ | ||
} | ||
|
||
Real | ||
BiharmonicLapBC::computeQpResidual() | ||
{ | ||
return -_lap_u.value(_t, _q_point[_qp]) * (_grad_test[_i][_qp] * _normals[_qp]); | ||
} |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,49 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
#include "FunctionPenaltyFluxBC.h" | ||
#include "Function.h" | ||
|
||
template <> | ||
InputParameters | ||
validParams<FunctionPenaltyFluxBC>() | ||
{ | ||
InputParameters params = validParams<IntegratedBC>(); | ||
params.addRequiredParam<Real>("penalty", "Penalty scalar"); | ||
params.addRequiredParam<FunctionName>("function", | ||
"Function used to compute the desired normal flux"); | ||
return params; | ||
} | ||
|
||
FunctionPenaltyFluxBC::FunctionPenaltyFluxBC(const InputParameters & parameters) | ||
: IntegratedBC(parameters), _func(getFunction("function")), _p(getParam<Real>("penalty")) | ||
{ | ||
} | ||
|
||
Real | ||
FunctionPenaltyFluxBC::computeQpResidual() | ||
{ | ||
Real dudn = _grad_u[_qp] * _normals[_qp]; | ||
Real dgdn = _func.gradient(_t, _q_point[_qp]) * _normals[_qp]; | ||
Real dvdn = _grad_test[_i][_qp] * _normals[_qp]; | ||
return _p * (dudn - dgdn) * dvdn; | ||
} | ||
|
||
Real | ||
FunctionPenaltyFluxBC::computeQpJacobian() | ||
{ | ||
Real dphi_j_dn = _grad_phi[_j][_qp] * _normals[_qp]; | ||
Real dphi_i_dn = _grad_test[_i][_qp] * _normals[_qp]; | ||
|
||
return _p * dphi_j_dn * dphi_i_dn; | ||
} |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,39 @@ | ||
/****************************************************************/ | ||
/* DO NOT MODIFY THIS HEADER */ | ||
/* MOOSE - Multiphysics Object Oriented Simulation Environment */ | ||
/* */ | ||
/* (c) 2010 Battelle Energy Alliance, LLC */ | ||
/* ALL RIGHTS RESERVED */ | ||
/* */ | ||
/* Prepared by Battelle Energy Alliance, LLC */ | ||
/* Under Contract No. DE-AC07-05ID14517 */ | ||
/* With the U. S. Department of Energy */ | ||
/* */ | ||
/* See COPYRIGHT for full restrictions */ | ||
/****************************************************************/ | ||
#include "Biharmonic.h" | ||
|
||
template <> | ||
InputParameters | ||
validParams<Biharmonic>() | ||
{ | ||
InputParameters params = validParams<Kernel>(); | ||
return params; | ||
} | ||
|
||
Biharmonic::Biharmonic(const InputParameters & parameters) | ||
: Kernel(parameters), _second_u(second()), _second_phi(secondPhi()), _second_test(secondTest()) | ||
{ | ||
} | ||
|
||
Real | ||
Biharmonic::computeQpResidual() | ||
{ | ||
return _second_u[_qp].tr() * _second_test[_i][_qp].tr(); | ||
} | ||
|
||
Real | ||
Biharmonic::computeQpJacobian() | ||
{ | ||
return _second_phi[_j][_qp].tr() * _second_test[_i][_qp].tr(); | ||
} |