LBPM/tests/TestInterfaceSpeed.cpp

152 lines
5.0 KiB
C++
Raw Normal View History

2013-12-11 19:50:21 -06:00
#include <iostream>
#include <math.h>
2015-02-18 13:41:23 -06:00
#include "analysis/TwoPhase.h"
2021-01-05 17:43:44 -06:00
#include "common/MPI.h"
#include "common/Communication.h"
2015-02-18 13:41:23 -06:00
#include "IO/Mesh.h"
#include "IO/Writer.h"
#include "ProfilerApp.h"
2013-12-11 19:50:21 -06:00
#define RADIUS 15
#define CAPRAD 20
#define HEIGHT 15.5
#define N 60
2014-01-27 10:43:24 -06:00
#define SPEED -1
2013-12-11 19:50:21 -06:00
#define PI 3.14159
int main (int argc, char *argv[])
{
2015-02-18 13:41:23 -06:00
// Initialize MPI
2020-10-08 10:03:42 -05:00
Utilities::startup( argc, argv );
2020-01-28 07:51:32 -06:00
Utilities::MPI comm( MPI_COMM_WORLD );
int rank = comm.getRank();
2021-01-06 10:58:43 -06:00
int toReturn = 0;
{
2018-05-15 15:03:29 -05:00
// Load inputs
2018-05-15 14:00:26 -05:00
string FILENAME = argv[1];
2018-05-15 14:06:28 -05:00
// Load inputs
if (rank==0) printf("Loading input database \n");
2018-05-15 15:19:11 -05:00
auto db = std::make_shared<Database>( FILENAME );
auto domain_db = db->getDatabase( "Domain" );
2018-05-15 14:06:28 -05:00
int Nx = domain_db->getVector<int>( "n" )[0];
int Ny = domain_db->getVector<int>( "n" )[1];
int Nz = domain_db->getVector<int>( "n" )[2];
2015-02-18 13:41:23 -06:00
2020-01-22 11:01:29 -06:00
auto Dm = std::make_shared<Domain>(domain_db,comm);
2018-05-15 15:11:35 -05:00
2018-05-24 20:10:08 -05:00
Nx+=2; Ny+=2; Nz+=2;
2021-01-05 17:43:44 -06:00
for (int i=0; i<Nx*Ny*Nz; i++) Dm->id[i] = 1;
2015-02-18 13:41:23 -06:00
2018-05-19 06:49:32 -05:00
Dm->CommInit();
2015-02-18 13:41:23 -06:00
2020-01-22 11:01:29 -06:00
auto Averages = std::make_shared<TwoPhase>(Dm);
2015-02-18 13:41:23 -06:00
int timestep=0;
2013-12-11 19:50:21 -06:00
double Cx,Cy,Cz;
double dist1,dist2;
2015-02-18 13:41:23 -06:00
2013-12-11 19:50:21 -06:00
Cx = Cy = Cz = N*0.5;
2021-01-05 17:43:44 -06:00
for (int k=0; k<Nz; k++){
for (int j=0; j<Ny; j++){
for (int i=0; i<Nx; i++){
2013-12-11 19:50:21 -06:00
dist2 = sqrt((i-Cx)*(i-Cx)+(j-Cy)*(j-Cy)+(k-Cz)*(k-Cz)) - CAPRAD;
dist2 = fabs(Cz-k)-HEIGHT;
2018-05-17 20:09:54 -05:00
Averages->Phase_tminus(i,j,k) = dist2;
2013-12-11 19:50:21 -06:00
}
}
}
Cz += SPEED;
2021-01-05 17:43:44 -06:00
for (int k=0; k<Nz; k++){
for (int j=0; j<Ny; j++){
for (int i=0; i<Nx; i++){
2013-12-11 19:50:21 -06:00
dist1 = sqrt((i-Cx)*(i-Cx)+(j-Cy)*(j-Cy)) - RADIUS;
dist2 = sqrt((i-Cx)*(i-Cx)+(j-Cy)*(j-Cy)+(k-Cz)*(k-Cz)) - CAPRAD;
dist2 = fabs(Cz-k)-HEIGHT;
2018-05-17 20:09:54 -05:00
Averages->SDs(i,j,k) = -dist1;
Averages->Phase(i,j,k) = dist2;
Averages->SDn(i,j,k) = dist2;
2013-12-11 19:50:21 -06:00
}
}
}
Cz += SPEED;
2021-01-05 17:43:44 -06:00
for (int k=0; k<Nz; k++){
for (int j=0; j<Ny; j++){
for (int i=0; i<Nx; i++){
2013-12-11 19:50:21 -06:00
dist2 = sqrt((i-Cx)*(i-Cx)+(j-Cy)*(j-Cy)+(k-Cz)*(k-Cz)) - CAPRAD;
dist2 = fabs(Cz-k)-HEIGHT;
2018-05-17 20:09:54 -05:00
Averages->Phase_tplus(i,j,k) = dist2;
2013-12-11 19:50:21 -06:00
}
}
}
//...........................................................................
2015-02-18 13:41:23 -06:00
//....................................................................
// The following only need to be done once
2018-05-17 20:09:54 -05:00
//Averages->SetupCubes(Dm);
Averages->UpdateSolid(); // unless the solid is deformable!
2015-02-18 13:41:23 -06:00
//....................................................................
// The following need to be called each time new averages are computed
2018-05-17 20:09:54 -05:00
Averages->Initialize();
Averages->UpdateMeshValues();
Averages->ComputeLocal();
Averages->Reduce();
Averages->PrintAll(timestep);
2015-02-18 13:41:23 -06:00
//....................................................................
2013-12-16 11:04:24 -06:00
2013-12-11 19:50:21 -06:00
printf("-------------------------------- \n");
2018-05-17 20:09:54 -05:00
printf("NWP volume = %f \n", Averages->nwp_volume);
printf("Area wn = %f, Analytical = %f \n", Averages->awn,2*PI*RADIUS*RADIUS);
printf("Area ns = %f, Analytical = %f \n", Averages->ans, 2*PI*RADIUS*(Nz-2)-4*PI*RADIUS*HEIGHT);
printf("Area ws = %f, Analytical = %f \n", Averages->aws, 4*PI*RADIUS*HEIGHT);
printf("Area s = %f, Analytical = %f \n", Averages->As, 2*PI*RADIUS*(Nz-2));
printf("Length wns = %f, Analytical = %f \n", Averages->lwns, 4*PI*RADIUS);
2021-12-14 15:46:00 -06:00
printf("Geodesic curvature (wn) = %f, Analytical = %f \n", Averages->KGwns_global, 0.0);
printf("Geodesic curvature (ws) = %f, Analytical = %f \n", Averages->KNwns_global, 4*PI);
2013-12-16 11:04:24 -06:00
// printf("Cos(theta_wns) = %f, Analytical = %f \n",efawns/lwns,1.0*RADIUS/CAPRAD);
2018-05-17 20:09:54 -05:00
printf("Interface Velocity = %f,%f,%f \n",Averages->vawn_global(0),Averages->vawn_global(1),Averages->vawn_global(2));
printf("Common Curve Velocity = %f,%f,%f \n",Averages->vawns_global(0),Averages->vawns_global(1),Averages->vawns_global(2));
2013-12-11 19:50:21 -06:00
printf("-------------------------------- \n");
//.........................................................................
2018-05-17 20:09:54 -05:00
if (fabs(Averages->awn - 2*PI*RADIUS*RADIUS)/(2*PI*RADIUS*RADIUS) > 0.02){
toReturn = 1;
printf("TestCylinderArea.cpp: error tolerance exceeded for wn area \n");
}
2018-05-17 20:09:54 -05:00
if (fabs(Averages->ans - (2*PI*RADIUS*(Nz-2)-4*PI*RADIUS*HEIGHT))/(2*PI*RADIUS*(Nz-2)-4*PI*RADIUS*HEIGHT)> 0.02 ){
toReturn = 2;
printf("TestCylinderArea.cpp: error tolerance exceeded for ns area \n");
}
2018-05-17 20:09:54 -05:00
if (fabs(Averages->aws - 4*PI*RADIUS*HEIGHT)/(4*PI*RADIUS*HEIGHT) > 0.02 ){
toReturn = 3;
printf("TestCylinderArea.cpp: error tolerance exceeded for ws area \n");
}
2018-05-17 20:09:54 -05:00
if (fabs(Averages->As - 2*PI*RADIUS*(Nz-2))/(2*PI*RADIUS*(Nz-2)) > 0.02 ){
toReturn = 4;
printf("TestCylinderArea.cpp: error tolerance exceeded for solid area \n");
}
2018-05-17 20:09:54 -05:00
if (fabs(Averages->lwns - 4*PI*RADIUS)/(4*PI*RADIUS) > 0.02 ){
toReturn = 5;
printf("TestCylinderArea.cpp: error tolerance exceeded for common curve length \n");
}
2018-05-24 20:10:08 -05:00
if ( fabs(Averages->vawn_global(2)+0.25) > 0.01){
2015-02-18 13:46:15 -06:00
printf("TestInterfaceSpeed: Error too high for kinematic velocity of wn interface \n");
toReturn = 6;
2015-02-18 13:46:15 -06:00
}
2018-05-24 20:10:08 -05:00
if ( fabs(Averages->vawns_global(2)+0.25) > 0.01){
2015-02-18 13:46:15 -06:00
printf("TestInterfaceSpeed: Error too high for kinematic velocity of common curve \n");
toReturn = 7;
2015-02-18 13:46:15 -06:00
}
2020-01-28 07:51:32 -06:00
comm.barrier();
2021-01-06 10:58:43 -06:00
}
Utilities::shutdown();
return toReturn;
2013-12-11 19:50:21 -06:00
}